gsl_randist__gauss.c 4.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144
  1. /* randist/gauss.c
  2. *
  3. * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2006, 2007 James Theiler, Brian Gough
  4. * Copyright (C) 2006 Charles Karney
  5. *
  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.
  10. *
  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.
  15. *
  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.
  19. */
  20. #include "gsl__config.h"
  21. #include <math.h>
  22. #include "gsl_math.h"
  23. #include "gsl_rng.h"
  24. #include "gsl_randist.h"
  25. /* Of the two methods provided below, I think the Polar method is more
  26. * efficient, but only when you are actually producing two random
  27. * deviates. We don't produce two, because then we'd have to save one
  28. * in a static variable for the next call, and that would screws up
  29. * re-entrant or threaded code, so we only produce one. This makes
  30. * the Ratio method suddenly more appealing.
  31. *
  32. * [Added by Charles Karney] We use Leva's implementation of the Ratio
  33. * method which avoids calling log() nearly all the time and makes the
  34. * Ratio method faster than the Polar method (when it produces just one
  35. * result per call). Timing per call (gcc -O2 on 866MHz Pentium,
  36. * average over 10^8 calls)
  37. *
  38. * Polar method: 660 ns
  39. * Ratio method: 368 ns
  40. *
  41. */
  42. /* Polar (Box-Mueller) method; See Knuth v2, 3rd ed, p122 */
  43. double
  44. gsl_ran_gaussian (const gsl_rng * r, const double sigma)
  45. {
  46. double x, y, r2;
  47. do
  48. {
  49. /* choose x,y in uniform square (-1,-1) to (+1,+1) */
  50. x = -1 + 2 * gsl_rng_uniform_pos (r);
  51. y = -1 + 2 * gsl_rng_uniform_pos (r);
  52. /* see if it is in the unit circle */
  53. r2 = x * x + y * y;
  54. }
  55. while (r2 > 1.0 || r2 == 0);
  56. /* Box-Muller transform */
  57. return sigma * y * sqrt (-2.0 * log (r2) / r2);
  58. }
  59. /* Ratio method (Kinderman-Monahan); see Knuth v2, 3rd ed, p130.
  60. * K+M, ACM Trans Math Software 3 (1977) 257-260.
  61. *
  62. * [Added by Charles Karney] This is an implementation of Leva's
  63. * modifications to the original K+M method; see:
  64. * J. L. Leva, ACM Trans Math Software 18 (1992) 449-453 and 454-455. */
  65. double
  66. gsl_ran_gaussian_ratio_method (const gsl_rng * r, const double sigma)
  67. {
  68. double u, v, x, y, Q;
  69. const double s = 0.449871; /* Constants from Leva */
  70. const double t = -0.386595;
  71. const double a = 0.19600;
  72. const double b = 0.25472;
  73. const double r1 = 0.27597;
  74. const double r2 = 0.27846;
  75. do /* This loop is executed 1.369 times on average */
  76. {
  77. /* Generate a point P = (u, v) uniform in a rectangle enclosing
  78. the K+M region v^2 <= - 4 u^2 log(u). */
  79. /* u in (0, 1] to avoid singularity at u = 0 */
  80. u = 1 - gsl_rng_uniform (r);
  81. /* v is in the asymmetric interval [-0.5, 0.5). However v = -0.5
  82. is rejected in the last part of the while clause. The
  83. resulting normal deviate is strictly symmetric about 0
  84. (provided that v is symmetric once v = -0.5 is excluded). */
  85. v = gsl_rng_uniform (r) - 0.5;
  86. /* Constant 1.7156 > sqrt(8/e) (for accuracy); but not by too
  87. much (for efficiency). */
  88. v *= 1.7156;
  89. /* Compute Leva's quadratic form Q */
  90. x = u - s;
  91. y = fabs (v) - t;
  92. Q = x * x + y * (a * y - b * x);
  93. /* Accept P if Q < r1 (Leva) */
  94. /* Reject P if Q > r2 (Leva) */
  95. /* Accept if v^2 <= -4 u^2 log(u) (K+M) */
  96. /* This final test is executed 0.012 times on average. */
  97. }
  98. while (Q >= r1 && (Q > r2 || v * v > -4 * u * u * log (u)));
  99. return sigma * (v / u); /* Return slope */
  100. }
  101. double
  102. gsl_ran_gaussian_pdf (const double x, const double sigma)
  103. {
  104. double u = x / fabs (sigma);
  105. double p = (1 / (sqrt (2 * M_PI) * fabs (sigma))) * exp (-u * u / 2);
  106. return p;
  107. }
  108. double
  109. gsl_ran_ugaussian (const gsl_rng * r)
  110. {
  111. return gsl_ran_gaussian (r, 1.0);
  112. }
  113. double
  114. gsl_ran_ugaussian_ratio_method (const gsl_rng * r)
  115. {
  116. return gsl_ran_gaussian_ratio_method (r, 1.0);
  117. }
  118. double
  119. gsl_ran_ugaussian_pdf (const double x)
  120. {
  121. return gsl_ran_gaussian_pdf (x, 1.0);
  122. }