3 * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2006, 2007 James Theiler, Brian Gough
4 * Copyright (C) 2006 Charles Karney
6 * This program is free software; you can redistribute it and/or modify
7 * it under the terms of the GNU General Public License as published by
8 * the Free Software Foundation; either version 3 of the License, or (at
9 * your option) any later version.
11 * This program is distributed in the hope that it will be useful, but
12 * WITHOUT ANY WARRANTY; without even the implied warranty of
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 * General Public License for more details.
16 * You should have received a copy of the GNU General Public License
17 * along with this program; if not, write to the Free Software
18 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
23 #include <gsl/gsl_math.h>
24 #include <gsl/gsl_rng.h>
25 #include <gsl/gsl_randist.h>
27 /* Of the two methods provided below, I think the Polar method is more
28 * efficient, but only when you are actually producing two random
29 * deviates. We don't produce two, because then we'd have to save one
30 * in a static variable for the next call, and that would screws up
31 * re-entrant or threaded code, so we only produce one. This makes
32 * the Ratio method suddenly more appealing.
34 * [Added by Charles Karney] We use Leva's implementation of the Ratio
35 * method which avoids calling log() nearly all the time and makes the
36 * Ratio method faster than the Polar method (when it produces just one
37 * result per call). Timing per call (gcc -O2 on 866MHz Pentium,
38 * average over 10^8 calls)
40 * Polar method: 660 ns
41 * Ratio method: 368 ns
45 /* Polar (Box-Mueller) method; See Knuth v2, 3rd ed, p122 */
48 gsl_ran_gaussian (const gsl_rng * r, const double sigma)
54 /* choose x,y in uniform square (-1,-1) to (+1,+1) */
55 x = -1 + 2 * gsl_rng_uniform_pos (r);
56 y = -1 + 2 * gsl_rng_uniform_pos (r);
58 /* see if it is in the unit circle */
61 while (r2 > 1.0 || r2 == 0);
63 /* Box-Muller transform */
64 return sigma * y * sqrt (-2.0 * log (r2) / r2);
67 /* Ratio method (Kinderman-Monahan); see Knuth v2, 3rd ed, p130.
68 * K+M, ACM Trans Math Software 3 (1977) 257-260.
70 * [Added by Charles Karney] This is an implementation of Leva's
71 * modifications to the original K+M method; see:
72 * J. L. Leva, ACM Trans Math Software 18 (1992) 449-453 and 454-455. */
75 gsl_ran_gaussian_ratio_method (const gsl_rng * r, const double sigma)
78 const double s = 0.449871; /* Constants from Leva */
79 const double t = -0.386595;
80 const double a = 0.19600;
81 const double b = 0.25472;
82 const double r1 = 0.27597;
83 const double r2 = 0.27846;
85 do /* This loop is executed 1.369 times on average */
87 /* Generate a point P = (u, v) uniform in a rectangle enclosing
88 the K+M region v^2 <= - 4 u^2 log(u). */
90 /* u in (0, 1] to avoid singularity at u = 0 */
91 u = 1 - gsl_rng_uniform (r);
93 /* v is in the asymmetric interval [-0.5, 0.5). However v = -0.5
94 is rejected in the last part of the while clause. The
95 resulting normal deviate is strictly symmetric about 0
96 (provided that v is symmetric once v = -0.5 is excluded). */
97 v = gsl_rng_uniform (r) - 0.5;
99 /* Constant 1.7156 > sqrt(8/e) (for accuracy); but not by too
100 much (for efficiency). */
103 /* Compute Leva's quadratic form Q */
106 Q = x * x + y * (a * y - b * x);
108 /* Accept P if Q < r1 (Leva) */
109 /* Reject P if Q > r2 (Leva) */
110 /* Accept if v^2 <= -4 u^2 log(u) (K+M) */
111 /* This final test is executed 0.012 times on average. */
113 while (Q >= r1 && (Q > r2 || v * v > -4 * u * u * log (u)));
115 return sigma * (v / u); /* Return slope */
119 gsl_ran_gaussian_pdf (const double x, const double sigma)
121 double u = x / fabs (sigma);
122 double p = (1 / (sqrt (2 * M_PI) * fabs (sigma))) * exp (-u * u / 2);
127 gsl_ran_ugaussian (const gsl_rng * r)
129 return gsl_ran_gaussian (r, 1.0);
133 gsl_ran_ugaussian_ratio_method (const gsl_rng * r)
135 return gsl_ran_gaussian_ratio_method (r, 1.0);
139 gsl_ran_ugaussian_pdf (const double x)
141 return gsl_ran_gaussian_pdf (x, 1.0);