3 * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2004, 2007 James Theiler, Brian Gough
5 * This program is free software; you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation; either version 3 of the License, or (at
8 * your option) any later version.
10 * This program is distributed in the hope that it will be useful, but
11 * WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13 * General Public License for more details.
15 * You should have received a copy of the GNU General Public License
16 * along with this program; if not, write to the Free Software
17 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
22 #include <gsl/gsl_rng.h>
23 #include <gsl/gsl_randist.h>
26 gsl_ran_dir_2d (const gsl_rng * r, double *x, double *y)
28 /* This method avoids trig, but it does take an average of 8/pi =
29 * 2.55 calls to the RNG, instead of one for the direct
30 * trigonometric method. */
35 u = -1 + 2 * gsl_rng_uniform (r);
36 v = -1 + 2 * gsl_rng_uniform (r);
39 while (s > 1.0 || s == 0.0);
41 /* This is the Von Neumann trick. See Knuth, v2, 3rd ed, p140
42 * (exercise 23). Note, no sin, cos, or sqrt ! */
44 *x = (u * u - v * v) / s;
47 /* Here is the more straightforward approach,
48 * s = sqrt (s); *x = u / s; *y = v / s;
49 * It has fewer total operations, but one of them is a sqrt */
53 gsl_ran_dir_2d_trig_method (const gsl_rng * r, double *x, double *y)
55 /* This is the obvious solution... */
56 /* It ain't clever, but since sin/cos are often hardware accelerated,
57 * it can be faster -- it is on my home Pentium -- than von Neumann's
58 * solution, or slower -- as it is on my Sun Sparc 20 at work
60 double t = 6.2831853071795864 * gsl_rng_uniform (r); /* 2*PI */
66 gsl_ran_dir_3d (const gsl_rng * r, double *x, double *y, double *z)
70 /* This is a variant of the algorithm for computing a random point
71 * on the unit sphere; the algorithm is suggested in Knuth, v2,
72 * 3rd ed, p136; and attributed to Robert E Knop, CACM, 13 (1970),
76 /* Begin with the polar method for getting x,y inside a unit circle
80 *x = -1 + 2 * gsl_rng_uniform (r);
81 *y = -1 + 2 * gsl_rng_uniform (r);
82 s = (*x) * (*x) + (*y) * (*y);
86 *z = -1 + 2 * s; /* z uniformly distributed from -1 to 1 */
87 a = 2 * sqrt (1 - s); /* factor to adjust x,y so that x^2+y^2
88 * is equal to 1-z^2 */
94 gsl_ran_dir_nd (const gsl_rng * r, size_t n, double *x)
98 /* See Knuth, v2, 3rd ed, p135-136. The method is attributed to
99 * G. W. Brown, in Modern Mathematics for the Engineer (1956).
100 * The idea is that gaussians G(x) have the property that
101 * G(x)G(y)G(z)G(...) is radially symmetric, a function only
102 * r = sqrt(x^2+y^2+...)
107 for (i = 0; i < n; ++i)
109 x[i] = gsl_ran_gaussian (r, 1.0);
115 for (i = 0; i < n; ++i)