gsl_randist__gamma.c 5.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221
  1. /* randist/gamma.c
  2. *
  3. * Copyright (C) 1996, 1997, 1998, 1999, 2000, 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_math.h"
  22. #include "gsl_sf_gamma.h"
  23. #include "gsl_rng.h"
  24. #include "gsl_randist.h"
  25. static double gamma_large (const gsl_rng * r, const double a);
  26. static double gamma_frac (const gsl_rng * r, const double a);
  27. /* The Gamma distribution of order a>0 is defined by:
  28. p(x) dx = {1 / \Gamma(a) b^a } x^{a-1} e^{-x/b} dx
  29. for x>0. If X and Y are independent gamma-distributed random
  30. variables of order a1 and a2 with the same scale parameter b, then
  31. X+Y has gamma distribution of order a1+a2.
  32. The algorithms below are from Knuth, vol 2, 2nd ed, p. 129. */
  33. double
  34. gsl_ran_gamma_knuth (const gsl_rng * r, const double a, const double b)
  35. {
  36. /* assume a > 0 */
  37. unsigned int na = floor (a);
  38. if (a == na)
  39. {
  40. return b * gsl_ran_gamma_int (r, na);
  41. }
  42. else if (na == 0)
  43. {
  44. return b * gamma_frac (r, a);
  45. }
  46. else
  47. {
  48. return b * (gsl_ran_gamma_int (r, na) + gamma_frac (r, a - na)) ;
  49. }
  50. }
  51. double
  52. gsl_ran_gamma_int (const gsl_rng * r, const unsigned int a)
  53. {
  54. if (a < 12)
  55. {
  56. unsigned int i;
  57. double prod = 1;
  58. for (i = 0; i < a; i++)
  59. {
  60. prod *= gsl_rng_uniform_pos (r);
  61. }
  62. /* Note: for 12 iterations we are safe against underflow, since
  63. the smallest positive random number is O(2^-32). This means
  64. the smallest possible product is 2^(-12*32) = 10^-116 which
  65. is within the range of double precision. */
  66. return -log (prod);
  67. }
  68. else
  69. {
  70. return gamma_large (r, (double) a);
  71. }
  72. }
  73. static double
  74. gamma_large (const gsl_rng * r, const double a)
  75. {
  76. /* Works only if a > 1, and is most efficient if a is large
  77. This algorithm, reported in Knuth, is attributed to Ahrens. A
  78. faster one, we are told, can be found in: J. H. Ahrens and
  79. U. Dieter, Computing 12 (1974) 223-246. */
  80. double sqa, x, y, v;
  81. sqa = sqrt (2 * a - 1);
  82. do
  83. {
  84. do
  85. {
  86. y = tan (M_PI * gsl_rng_uniform (r));
  87. x = sqa * y + a - 1;
  88. }
  89. while (x <= 0);
  90. v = gsl_rng_uniform (r);
  91. }
  92. while (v > (1 + y * y) * exp ((a - 1) * log (x / (a - 1)) - sqa * y));
  93. return x;
  94. }
  95. static double
  96. gamma_frac (const gsl_rng * r, const double a)
  97. {
  98. /* This is exercise 16 from Knuth; see page 135, and the solution is
  99. on page 551. */
  100. double p, q, x, u, v;
  101. p = M_E / (a + M_E);
  102. do
  103. {
  104. u = gsl_rng_uniform (r);
  105. v = gsl_rng_uniform_pos (r);
  106. if (u < p)
  107. {
  108. x = exp ((1 / a) * log (v));
  109. q = exp (-x);
  110. }
  111. else
  112. {
  113. x = 1 - log (v);
  114. q = exp ((a - 1) * log (x));
  115. }
  116. }
  117. while (gsl_rng_uniform (r) >= q);
  118. return x;
  119. }
  120. double
  121. gsl_ran_gamma_pdf (const double x, const double a, const double b)
  122. {
  123. if (x < 0)
  124. {
  125. return 0 ;
  126. }
  127. else if (x == 0)
  128. {
  129. if (a == 1)
  130. return 1/b ;
  131. else
  132. return 0 ;
  133. }
  134. else if (a == 1)
  135. {
  136. return exp(-x/b)/b ;
  137. }
  138. else
  139. {
  140. double p;
  141. double lngamma = gsl_sf_lngamma (a);
  142. p = exp ((a - 1) * log (x/b) - x/b - lngamma)/b;
  143. return p;
  144. }
  145. }
  146. /* New version based on Marsaglia and Tsang, "A Simple Method for
  147. * generating gamma variables", ACM Transactions on Mathematical
  148. * Software, Vol 26, No 3 (2000), p363-372.
  149. *
  150. * Implemented by J.D.Lamb@btinternet.com, minor modifications for GSL
  151. * by Brian Gough
  152. */
  153. double
  154. gsl_ran_gamma_mt (const gsl_rng * r, const double a, const double b)
  155. {
  156. return gsl_ran_gamma (r, a, b);
  157. }
  158. double
  159. gsl_ran_gamma (const gsl_rng * r, const double a, const double b)
  160. {
  161. /* assume a > 0 */
  162. if (a < 1)
  163. {
  164. double u = gsl_rng_uniform_pos (r);
  165. return gsl_ran_gamma (r, 1.0 + a, b) * pow (u, 1.0 / a);
  166. }
  167. {
  168. double x, v, u;
  169. double d = a - 1.0 / 3.0;
  170. double c = (1.0 / 3.0) / sqrt (d);
  171. while (1)
  172. {
  173. do
  174. {
  175. x = gsl_ran_gaussian_ziggurat (r, 1.0);
  176. v = 1.0 + c * x;
  177. }
  178. while (v <= 0);
  179. v = v * v * v;
  180. u = gsl_rng_uniform_pos (r);
  181. if (u < 1 - 0.0331 * x * x * x * x)
  182. break;
  183. if (log (u) < 0.5 * x * x + d * (1 - v + log (v)))
  184. break;
  185. }
  186. return b * d * v;
  187. }
  188. }