gsl_randist__binomial_tpe.c 13 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382
  1. /* randist/binomial_tpe.c
  2. *
  3. * Copyright (C) 1996, 2003, 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. #include "gsl_pow_int.h"
  24. #include "gsl_sf_gamma.h"
  25. /* The binomial distribution has the form,
  26. f(x) = n!/(x!(n-x)!) * p^x (1-p)^(n-x) for integer 0 <= x <= n
  27. = 0 otherwise
  28. This implementation follows the public domain ranlib function
  29. "ignbin", the bulk of which is the BTPE (Binomial Triangle
  30. Parallelogram Exponential) algorithm introduced in
  31. Kachitvichyanukul and Schmeiser[1]. It has been translated to use
  32. modern C coding standards.
  33. If n is small and/or p is near 0 or near 1 (specifically, if
  34. n*min(p,1-p) < SMALL_MEAN), then a different algorithm, called
  35. BINV, is used which has an average runtime that scales linearly
  36. with n*min(p,1-p).
  37. But for larger problems, the BTPE algorithm takes the form of two
  38. functions b(x) and t(x) -- "bottom" and "top" -- for which b(x) <
  39. f(x)/f(M) < t(x), with M = floor(n*p+p). b(x) defines a triangular
  40. region, and t(x) includes a parallelogram and two tails. Details
  41. (including a nice drawing) are in the paper.
  42. [1] Kachitvichyanukul, V. and Schmeiser, B. W. Binomial Random
  43. Variate Generation. Communications of the ACM, 31, 2 (February,
  44. 1988) 216.
  45. Note, Bruce Schmeiser (personal communication) points out that if
  46. you want very fast binomial deviates, and you are happy with
  47. approximate results, and/or n and n*p are both large, then you can
  48. just use gaussian estimates: mean=n*p, variance=n*p*(1-p).
  49. This implementation by James Theiler, April 2003, after obtaining
  50. permission -- and some good advice -- from Drs. Kachitvichyanukul
  51. and Schmeiser to use their code as a starting point, and then doing
  52. a little bit of tweaking.
  53. Additional polishing for GSL coding standards by Brian Gough. */
  54. #define SMALL_MEAN 14 /* If n*p < SMALL_MEAN then use BINV
  55. algorithm. The ranlib
  56. implementation used cutoff=30; but
  57. on my computer 14 works better */
  58. #define BINV_CUTOFF 110 /* In BINV, do not permit ix too large */
  59. #define FAR_FROM_MEAN 20 /* If ix-n*p is larger than this, then
  60. use the "squeeze" algorithm.
  61. Ranlib used 20, and this seems to
  62. be the best choice on my machine as
  63. well */
  64. #define LNFACT(x) gsl_sf_lnfact(x)
  65. inline static double
  66. Stirling (double y1)
  67. {
  68. double y2 = y1 * y1;
  69. double s =
  70. (13860.0 -
  71. (462.0 - (132.0 - (99.0 - 140.0 / y2) / y2) / y2) / y2) / y1 / 166320.0;
  72. return s;
  73. }
  74. unsigned int
  75. gsl_ran_binomial_tpe (const gsl_rng * rng, double p, unsigned int n)
  76. {
  77. return gsl_ran_binomial (rng, p, n);
  78. }
  79. unsigned int
  80. gsl_ran_binomial (const gsl_rng * rng, double p, unsigned int n)
  81. {
  82. int ix; /* return value */
  83. int flipped = 0;
  84. double q, s, np;
  85. if (n == 0)
  86. return 0;
  87. if (p > 0.5)
  88. {
  89. p = 1.0 - p; /* work with small p */
  90. flipped = 1;
  91. }
  92. q = 1 - p;
  93. s = p / q;
  94. np = n * p;
  95. /* Inverse cdf logic for small mean (BINV in K+S) */
  96. if (np < SMALL_MEAN)
  97. {
  98. double f0 = gsl_pow_int (q, n); /* f(x), starting with x=0 */
  99. while (1)
  100. {
  101. /* This while(1) loop will almost certainly only loop once; but
  102. * if u=1 to within a few epsilons of machine precision, then it
  103. * is possible for roundoff to prevent the main loop over ix to
  104. * achieve its proper value. following the ranlib implementation,
  105. * we introduce a check for that situation, and when it occurs,
  106. * we just try again.
  107. */
  108. double f = f0;
  109. double u = gsl_rng_uniform (rng);
  110. for (ix = 0; ix <= BINV_CUTOFF; ++ix)
  111. {
  112. if (u < f)
  113. goto Finish;
  114. u -= f;
  115. /* Use recursion f(x+1) = f(x)*[(n-x)/(x+1)]*[p/(1-p)] */
  116. f *= s * (n - ix) / (ix + 1);
  117. }
  118. /* It should be the case that the 'goto Finish' was encountered
  119. * before this point was ever reached. But if we have reached
  120. * this point, then roundoff has prevented u from decreasing
  121. * all the way to zero. This can happen only if the initial u
  122. * was very nearly equal to 1, which is a rare situation. In
  123. * that rare situation, we just try again.
  124. *
  125. * Note, following the ranlib implementation, we loop ix only to
  126. * a hardcoded value of SMALL_MEAN_LARGE_N=110; we could have
  127. * looped to n, and 99.99...% of the time it won't matter. This
  128. * choice, I think is a little more robust against the rare
  129. * roundoff error. If n>LARGE_N, then it is technically
  130. * possible for ix>LARGE_N, but it is astronomically rare, and
  131. * if ix is that large, it is more likely due to roundoff than
  132. * probability, so better to nip it at LARGE_N than to take a
  133. * chance that roundoff will somehow conspire to produce an even
  134. * larger (and more improbable) ix. If n<LARGE_N, then once
  135. * ix=n, f=0, and the loop will continue until ix=LARGE_N.
  136. */
  137. }
  138. }
  139. else
  140. {
  141. /* For n >= SMALL_MEAN, we invoke the BTPE algorithm */
  142. int k;
  143. double ffm = np + p; /* ffm = n*p+p */
  144. int m = (int) ffm; /* m = int floor[n*p+p] */
  145. double fm = m; /* fm = double m; */
  146. double xm = fm + 0.5; /* xm = half integer mean (tip of triangle) */
  147. double npq = np * q; /* npq = n*p*q */
  148. /* Compute cumulative area of tri, para, exp tails */
  149. /* p1: radius of triangle region; since height=1, also: area of region */
  150. /* p2: p1 + area of parallelogram region */
  151. /* p3: p2 + area of left tail */
  152. /* p4: p3 + area of right tail */
  153. /* pi/p4: probability of i'th area (i=1,2,3,4) */
  154. /* Note: magic numbers 2.195, 4.6, 0.134, 20.5, 15.3 */
  155. /* These magic numbers are not adjustable...at least not easily! */
  156. double p1 = floor (2.195 * sqrt (npq) - 4.6 * q) + 0.5;
  157. /* xl, xr: left and right edges of triangle */
  158. double xl = xm - p1;
  159. double xr = xm + p1;
  160. /* Parameter of exponential tails */
  161. /* Left tail: t(x) = c*exp(-lambda_l*[xl - (x+0.5)]) */
  162. /* Right tail: t(x) = c*exp(-lambda_r*[(x+0.5) - xr]) */
  163. double c = 0.134 + 20.5 / (15.3 + fm);
  164. double p2 = p1 * (1.0 + c + c);
  165. double al = (ffm - xl) / (ffm - xl * p);
  166. double lambda_l = al * (1.0 + 0.5 * al);
  167. double ar = (xr - ffm) / (xr * q);
  168. double lambda_r = ar * (1.0 + 0.5 * ar);
  169. double p3 = p2 + c / lambda_l;
  170. double p4 = p3 + c / lambda_r;
  171. double var, accept;
  172. double u, v; /* random variates */
  173. TryAgain:
  174. /* generate random variates, u specifies which region: Tri, Par, Tail */
  175. u = gsl_rng_uniform (rng) * p4;
  176. v = gsl_rng_uniform (rng);
  177. if (u <= p1)
  178. {
  179. /* Triangular region */
  180. ix = (int) (xm - p1 * v + u);
  181. goto Finish;
  182. }
  183. else if (u <= p2)
  184. {
  185. /* Parallelogram region */
  186. double x = xl + (u - p1) / c;
  187. v = v * c + 1.0 - fabs (x - xm) / p1;
  188. if (v > 1.0 || v <= 0.0)
  189. goto TryAgain;
  190. ix = (int) x;
  191. }
  192. else if (u <= p3)
  193. {
  194. /* Left tail */
  195. ix = (int) (xl + log (v) / lambda_l);
  196. if (ix < 0)
  197. goto TryAgain;
  198. v *= ((u - p2) * lambda_l);
  199. }
  200. else
  201. {
  202. /* Right tail */
  203. ix = (int) (xr - log (v) / lambda_r);
  204. if (ix > (double) n)
  205. goto TryAgain;
  206. v *= ((u - p3) * lambda_r);
  207. }
  208. /* At this point, the goal is to test whether v <= f(x)/f(m)
  209. *
  210. * v <= f(x)/f(m) = (m!(n-m)! / (x!(n-x)!)) * (p/q)^{x-m}
  211. *
  212. */
  213. /* Here is a direct test using logarithms. It is a little
  214. * slower than the various "squeezing" computations below, but
  215. * if things are working, it should give exactly the same answer
  216. * (given the same random number seed). */
  217. #ifdef DIRECT
  218. var = log (v);
  219. accept =
  220. LNFACT (m) + LNFACT (n - m) - LNFACT (ix) - LNFACT (n - ix)
  221. + (ix - m) * log (p / q);
  222. #else /* SQUEEZE METHOD */
  223. /* More efficient determination of whether v < f(x)/f(M) */
  224. k = abs (ix - m);
  225. if (k <= FAR_FROM_MEAN)
  226. {
  227. /*
  228. * If ix near m (ie, |ix-m|<FAR_FROM_MEAN), then do
  229. * explicit evaluation using recursion relation for f(x)
  230. */
  231. double g = (n + 1) * s;
  232. double f = 1.0;
  233. var = v;
  234. if (m < ix)
  235. {
  236. int i;
  237. for (i = m + 1; i <= ix; i++)
  238. {
  239. f *= (g / i - s);
  240. }
  241. }
  242. else if (m > ix)
  243. {
  244. int i;
  245. for (i = ix + 1; i <= m; i++)
  246. {
  247. f /= (g / i - s);
  248. }
  249. }
  250. accept = f;
  251. }
  252. else
  253. {
  254. /* If ix is far from the mean m: k=ABS(ix-m) large */
  255. var = log (v);
  256. if (k < npq / 2 - 1)
  257. {
  258. /* "Squeeze" using upper and lower bounds on
  259. * log(f(x)) The squeeze condition was derived
  260. * under the condition k < npq/2-1 */
  261. double amaxp =
  262. k / npq * ((k * (k / 3.0 + 0.625) + (1.0 / 6.0)) / npq + 0.5);
  263. double ynorm = -(k * k / (2.0 * npq));
  264. if (var < ynorm - amaxp)
  265. goto Finish;
  266. if (var > ynorm + amaxp)
  267. goto TryAgain;
  268. }
  269. /* Now, again: do the test log(v) vs. log f(x)/f(M) */
  270. #if USE_EXACT
  271. /* This is equivalent to the above, but is a little (~20%) slower */
  272. /* There are five log's vs three above, maybe that's it? */
  273. accept = LNFACT (m) + LNFACT (n - m)
  274. - LNFACT (ix) - LNFACT (n - ix) + (ix - m) * log (p / q);
  275. #else /* USE STIRLING */
  276. /* The "#define Stirling" above corresponds to the first five
  277. * terms in asymptoic formula for
  278. * log Gamma (y) - (y-0.5)log(y) + y - 0.5 log(2*pi);
  279. * See Abramowitz and Stegun, eq 6.1.40
  280. */
  281. /* Note below: two Stirling's are added, and two are
  282. * subtracted. In both K+S, and in the ranlib
  283. * implementation, all four are added. I (jt) believe that
  284. * is a mistake -- this has been confirmed by personal
  285. * correspondence w/ Dr. Kachitvichyanukul. Note, however,
  286. * the corrections are so small, that I couldn't find an
  287. * example where it made a difference that could be
  288. * observed, let alone tested. In fact, define'ing Stirling
  289. * to be zero gave identical results!! In practice, alv is
  290. * O(1), ranging 0 to -10 or so, while the Stirling
  291. * correction is typically O(10^{-5}) ...setting the
  292. * correction to zero gives about a 2% performance boost;
  293. * might as well keep it just to be pendantic. */
  294. {
  295. double x1 = ix + 1.0;
  296. double w1 = n - ix + 1.0;
  297. double f1 = fm + 1.0;
  298. double z1 = n + 1.0 - fm;
  299. accept = xm * log (f1 / x1) + (n - m + 0.5) * log (z1 / w1)
  300. + (ix - m) * log (w1 * p / (x1 * q))
  301. + Stirling (f1) + Stirling (z1) - Stirling (x1) - Stirling (w1);
  302. }
  303. #endif
  304. #endif
  305. }
  306. if (var <= accept)
  307. {
  308. goto Finish;
  309. }
  310. else
  311. {
  312. goto TryAgain;
  313. }
  314. }
  315. Finish:
  316. return (flipped) ? (n - ix) : (unsigned int)ix;
  317. }