gsl_randist__sphere.c 3.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120
  1. /* randist/sphere.c
  2. *
  3. * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2004, 2007 James Theiler, Brian Gough
  4. *
  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.
  9. *
  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.
  14. *
  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.
  18. */
  19. #include "gsl__config.h"
  20. #include <math.h>
  21. #include "gsl_rng.h"
  22. #include "gsl_randist.h"
  23. void
  24. gsl_ran_dir_2d (const gsl_rng * r, double *x, double *y)
  25. {
  26. /* This method avoids trig, but it does take an average of 8/pi =
  27. * 2.55 calls to the RNG, instead of one for the direct
  28. * trigonometric method. */
  29. double u, v, s;
  30. do
  31. {
  32. u = -1 + 2 * gsl_rng_uniform (r);
  33. v = -1 + 2 * gsl_rng_uniform (r);
  34. s = u * u + v * v;
  35. }
  36. while (s > 1.0 || s == 0.0);
  37. /* This is the Von Neumann trick. See Knuth, v2, 3rd ed, p140
  38. * (exercise 23). Note, no sin, cos, or sqrt ! */
  39. *x = (u * u - v * v) / s;
  40. *y = 2 * u * v / s;
  41. /* Here is the more straightforward approach,
  42. * s = sqrt (s); *x = u / s; *y = v / s;
  43. * It has fewer total operations, but one of them is a sqrt */
  44. }
  45. void
  46. gsl_ran_dir_2d_trig_method (const gsl_rng * r, double *x, double *y)
  47. {
  48. /* This is the obvious solution... */
  49. /* It ain't clever, but since sin/cos are often hardware accelerated,
  50. * it can be faster -- it is on my home Pentium -- than von Neumann's
  51. * solution, or slower -- as it is on my Sun Sparc 20 at work
  52. */
  53. double t = 6.2831853071795864 * gsl_rng_uniform (r); /* 2*PI */
  54. *x = cos (t);
  55. *y = sin (t);
  56. }
  57. void
  58. gsl_ran_dir_3d (const gsl_rng * r, double *x, double *y, double *z)
  59. {
  60. double s, a;
  61. /* This is a variant of the algorithm for computing a random point
  62. * on the unit sphere; the algorithm is suggested in Knuth, v2,
  63. * 3rd ed, p136; and attributed to Robert E Knop, CACM, 13 (1970),
  64. * 326.
  65. */
  66. /* Begin with the polar method for getting x,y inside a unit circle
  67. */
  68. do
  69. {
  70. *x = -1 + 2 * gsl_rng_uniform (r);
  71. *y = -1 + 2 * gsl_rng_uniform (r);
  72. s = (*x) * (*x) + (*y) * (*y);
  73. }
  74. while (s > 1.0);
  75. *z = -1 + 2 * s; /* z uniformly distributed from -1 to 1 */
  76. a = 2 * sqrt (1 - s); /* factor to adjust x,y so that x^2+y^2
  77. * is equal to 1-z^2 */
  78. *x *= a;
  79. *y *= a;
  80. }
  81. void
  82. gsl_ran_dir_nd (const gsl_rng * r, size_t n, double *x)
  83. {
  84. double d;
  85. size_t i;
  86. /* See Knuth, v2, 3rd ed, p135-136. The method is attributed to
  87. * G. W. Brown, in Modern Mathematics for the Engineer (1956).
  88. * The idea is that gaussians G(x) have the property that
  89. * G(x)G(y)G(z)G(...) is radially symmetric, a function only
  90. * r = sqrt(x^2+y^2+...)
  91. */
  92. d = 0;
  93. do
  94. {
  95. for (i = 0; i < n; ++i)
  96. {
  97. x[i] = gsl_ran_gaussian (r, 1.0);
  98. d += x[i] * x[i];
  99. }
  100. }
  101. while (d == 0);
  102. d = sqrt (d);
  103. for (i = 0; i < n; ++i)
  104. {
  105. x[i] /= d;
  106. }
  107. }