gsl_randist__discrete.c 13 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399
  1. /* randist/discrete.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. /*
  20. Random Discrete Events
  21. Given K discrete events with different probabilities P[k]
  22. produce a value k consistent with its probability.
  23. This program is free software; you can redistribute it and/or
  24. modify it under the terms of the GNU General Public License as
  25. published by the Free Software Foundation; either version 3 of the
  26. License, or (at your option) any later version.
  27. This program is distributed in the hope that it will be useful, but
  28. WITHOUT ANY WARRANTY; without even the implied warranty of
  29. MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
  30. General Public License for more details. You should have received
  31. a copy of the GNU General Public License along with this program;
  32. if not, write to the Free Foundation, Inc., 59 Temple Place, Suite
  33. 330, Boston, MA 02111-1307 USA
  34. */
  35. /*
  36. * Based on: Alastair J Walker, An efficient method for generating
  37. * discrete random variables with general distributions, ACM Trans
  38. * Math Soft 3, 253-256 (1977). See also: D. E. Knuth, The Art of
  39. * Computer Programming, Volume 2 (Seminumerical algorithms), 3rd
  40. * edition, Addison-Wesley (1997), p120.
  41. * Walker's algorithm does some preprocessing, and provides two
  42. * arrays: floating point F[k] and integer A[k]. A value k is chosen
  43. * from 0..K-1 with equal likelihood, and then a uniform random number
  44. * u is compared to F[k]. If it is less than F[k], then k is
  45. * returned. Otherwise, A[k] is returned.
  46. * Walker's original paper describes an O(K^2) algorithm for setting
  47. * up the F and A arrays. I found this disturbing since I wanted to
  48. * use very large values of K. I'm sure I'm not the first to realize
  49. * this, but in fact the preprocessing can be done in O(K) steps.
  50. * A figure of merit for the preprocessing is the average value for
  51. * the F[k]'s (that is, SUM_k F[k]/K); this corresponds to the
  52. * probability that k is returned, instead of A[k], thereby saving a
  53. * redirection. Walker's O(K^2) preprocessing will generally improve
  54. * that figure of merit, compared to my cheaper O(K) method; from some
  55. * experiments with a perl script, I get values of around 0.6 for my
  56. * method and just under 0.75 for Walker's. Knuth has pointed out
  57. * that finding _the_ optimum lookup tables, which maximize the
  58. * average F[k], is a combinatorially difficult problem. But any
  59. * valid preprocessing will still provide O(1) time for the call to
  60. * gsl_ran_discrete(). I find that if I artificially set F[k]=1 --
  61. * ie, better than optimum! -- I get a speedup of maybe 20%, so that's
  62. * the maximum I could expect from the most expensive preprocessing.
  63. * Folding in the difference of 0.6 vs 0.75, I'd estimate that the
  64. * speedup would be less than 10%.
  65. * I've not implemented it here, but one compromise is to sort the
  66. * probabilities once, and then work from the two ends inward. This
  67. * requires O(K log K), still lots cheaper than O(K^2), and from my
  68. * experiments with the perl script, the figure of merit is within
  69. * about 0.01 for K up to 1000, and no sign of diverging (in fact,
  70. * they seemed to be converging, but it's hard to say with just a
  71. * handful of runs).
  72. * The O(K) algorithm goes through all the p_k's and decides if they
  73. * are "smalls" or "bigs" according to whether they are less than or
  74. * greater than the mean value 1/K. The indices to the smalls and the
  75. * bigs are put in separate stacks, and then we work through the
  76. * stacks together. For each small, we pair it up with the next big
  77. * in the stack (Walker always wanted to pair up the smallest small
  78. * with the biggest big). The small "borrows" from the big just
  79. * enough to bring the small up to mean. This reduces the size of the
  80. * big, so the (smaller) big is compared again to the mean, and if it
  81. * is smaller, it gets "popped" from the big stack and "pushed" to the
  82. * small stack. Otherwise, it stays put. Since every time we pop a
  83. * small, we are able to deal with it right then and there, and we
  84. * never have to pop more than K smalls, then the algorithm is O(K).
  85. * This implementation sets up two separate stacks, and allocates K
  86. * elements between them. Since neither stack ever grows, we do an
  87. * extra O(K) pass through the data to determine how many smalls and
  88. * bigs there are to begin with and allocate appropriately. In all
  89. * there are 2*K*sizeof(double) transient bytes of memory that are
  90. * used than returned, and K*(sizeof(int)+sizeof(double)) bytes used
  91. * in the lookup table.
  92. * Walker spoke of using two random numbers (an integer 0..K-1, and a
  93. * floating point u in [0,1]), but Knuth points out that one can just
  94. * use the integer and fractional parts of K*u where u is in [0,1].
  95. * In fact, Knuth further notes that taking F'[k]=(k+F[k])/K, one can
  96. * directly compare u to F'[k] without having to explicitly set
  97. * u=K*u-int(K*u).
  98. * Usage:
  99. * Starting with an array of probabilities P, initialize and do
  100. * preprocessing with a call to:
  101. * gsl_rng *r;
  102. * gsl_ran_discrete_t *f;
  103. * f = gsl_ran_discrete_preproc(K,P);
  104. * Then, whenever a random index 0..K-1 is desired, use
  105. * k = gsl_ran_discrete(r,f);
  106. * Note that several different randevent struct's can be
  107. * simultaneously active.
  108. * Aside: A very clever alternative approach is described in
  109. * Abramowitz and Stegun, p 950, citing: Marsaglia, Random variables
  110. * and computers, Proc Third Prague Conference in Probability Theory,
  111. * 1962. A more accesible reference is: G. Marsaglia, Generating
  112. * discrete random numbers in a computer, Comm ACM 6, 37-38 (1963).
  113. * If anybody is interested, I (jt) have also coded up this version as
  114. * part of another software package. However, I've done some
  115. * comparisons, and the Walker method is both faster and more stingy
  116. * with memory. So, in the end I decided not to include it with the
  117. * GSL package.
  118. * Written 26 Jan 1999, James Theiler, jt@lanl.gov
  119. * Adapted to GSL, 30 Jan 1999, jt
  120. */
  121. #include "gsl__config.h"
  122. #include <stdio.h>
  123. #include <stdlib.h>
  124. #include <math.h>
  125. #include "gsl_errno.h"
  126. #include "gsl_rng.h"
  127. #include "gsl_randist.h"
  128. #define DEBUG 0
  129. #define KNUTH_CONVENTION 1 /* Saves a few steps of arithmetic
  130. * in the call to gsl_ran_discrete()
  131. */
  132. /*** Begin Stack (this code is used just in this file) ***/
  133. /* Stack code converted to use unsigned indices (i.e. s->i == 0 now
  134. means an empty stack, instead of -1), for consistency and to give a
  135. bigger allowable range. BJG */
  136. typedef struct {
  137. size_t N; /* max number of elts on stack */
  138. size_t *v; /* array of values on the stack */
  139. size_t i; /* index of top of stack */
  140. } gsl_stack_t;
  141. static gsl_stack_t *
  142. new_stack(size_t N) {
  143. gsl_stack_t *s;
  144. s = (gsl_stack_t *)malloc(sizeof(gsl_stack_t));
  145. s->N = N;
  146. s->i = 0; /* indicates stack is empty */
  147. s->v = (size_t *)malloc(sizeof(size_t)*N);
  148. return s;
  149. }
  150. static void
  151. push_stack(gsl_stack_t *s, size_t v)
  152. {
  153. if ((s->i) >= (s->N)) {
  154. fprintf(stderr,"Cannot push stack!\n");
  155. abort(); /* FIXME: fatal!! */
  156. }
  157. (s->v)[s->i] = v;
  158. s->i += 1;
  159. }
  160. static size_t pop_stack(gsl_stack_t *s)
  161. {
  162. if ((s->i) == 0) {
  163. fprintf(stderr,"Cannot pop stack!\n");
  164. abort(); /* FIXME: fatal!! */
  165. }
  166. s->i -= 1;
  167. return ((s->v)[s->i]);
  168. }
  169. static inline size_t size_stack(const gsl_stack_t *s)
  170. {
  171. return s->i;
  172. }
  173. static void free_stack(gsl_stack_t *s)
  174. {
  175. free((char *)(s->v));
  176. free((char *)s);
  177. }
  178. /*** End Stack ***/
  179. /*** Begin Walker's Algorithm ***/
  180. gsl_ran_discrete_t *
  181. gsl_ran_discrete_preproc(size_t Kevents, const double *ProbArray)
  182. {
  183. size_t k,b,s;
  184. gsl_ran_discrete_t *g;
  185. size_t nBigs, nSmalls;
  186. gsl_stack_t *Bigs;
  187. gsl_stack_t *Smalls;
  188. double *E;
  189. double pTotal = 0.0, mean, d;
  190. if (Kevents < 1) {
  191. /* Could probably treat Kevents=1 as a special case */
  192. GSL_ERROR_VAL ("number of events must be a positive integer",
  193. GSL_EINVAL, 0);
  194. }
  195. /* Make sure elements of ProbArray[] are positive.
  196. * Won't enforce that sum is unity; instead will just normalize
  197. */
  198. for (k=0; k<Kevents; ++k) {
  199. if (ProbArray[k] < 0) {
  200. GSL_ERROR_VAL ("probabilities must be non-negative",
  201. GSL_EINVAL, 0) ;
  202. }
  203. pTotal += ProbArray[k];
  204. }
  205. /* Begin setting up the main "object" (just a struct, no steroids) */
  206. g = (gsl_ran_discrete_t *)malloc(sizeof(gsl_ran_discrete_t));
  207. g->K = Kevents;
  208. g->F = (double *)malloc(sizeof(double)*Kevents);
  209. g->A = (size_t *)malloc(sizeof(size_t)*Kevents);
  210. E = (double *)malloc(sizeof(double)*Kevents);
  211. if (E==NULL) {
  212. GSL_ERROR_VAL ("Cannot allocate memory for randevent", GSL_ENOMEM, 0);
  213. }
  214. for (k=0; k<Kevents; ++k) {
  215. E[k] = ProbArray[k]/pTotal;
  216. }
  217. /* Now create the Bigs and the Smalls */
  218. mean = 1.0/Kevents;
  219. nSmalls=nBigs=0;
  220. for (k=0; k<Kevents; ++k) {
  221. if (E[k] < mean) ++nSmalls;
  222. else ++nBigs;
  223. }
  224. Bigs = new_stack(nBigs);
  225. Smalls = new_stack(nSmalls);
  226. for (k=0; k<Kevents; ++k) {
  227. if (E[k] < mean) {
  228. push_stack(Smalls,k);
  229. }
  230. else {
  231. push_stack(Bigs,k);
  232. }
  233. }
  234. /* Now work through the smalls */
  235. while (size_stack(Smalls) > 0) {
  236. s = pop_stack(Smalls);
  237. if (size_stack(Bigs) == 0) {
  238. (g->A)[s]=s;
  239. (g->F)[s]=1.0;
  240. continue;
  241. }
  242. b = pop_stack(Bigs);
  243. (g->A)[s]=b;
  244. (g->F)[s]=Kevents*E[s];
  245. #if DEBUG
  246. fprintf(stderr,"s=%2d, A=%2d, F=%.4f\n",s,(g->A)[s],(g->F)[s]);
  247. #endif
  248. d = mean - E[s];
  249. E[s] += d; /* now E[s] == mean */
  250. E[b] -= d;
  251. if (E[b] < mean) {
  252. push_stack(Smalls,b); /* no longer big, join ranks of the small */
  253. }
  254. else if (E[b] > mean) {
  255. push_stack(Bigs,b); /* still big, put it back where you found it */
  256. }
  257. else {
  258. /* E[b]==mean implies it is finished too */
  259. (g->A)[b]=b;
  260. (g->F)[b]=1.0;
  261. }
  262. }
  263. while (size_stack(Bigs) > 0) {
  264. b = pop_stack(Bigs);
  265. (g->A)[b]=b;
  266. (g->F)[b]=1.0;
  267. }
  268. /* Stacks have been emptied, and A and F have been filled */
  269. if ( size_stack(Smalls) != 0) {
  270. GSL_ERROR_VAL ("Smalls stack has not been emptied",
  271. GSL_ESANITY, 0 );
  272. }
  273. #if 0
  274. /* if 1, then artificially set all F[k]'s to unity. This will
  275. * give wrong answers, but you'll get them faster. But, not
  276. * that much faster (I get maybe 20%); that's an upper bound
  277. * on what the optimal preprocessing would give.
  278. */
  279. for (k=0; k<Kevents; ++k) {
  280. (g->F)[k] = 1.0;
  281. }
  282. #endif
  283. #if KNUTH_CONVENTION
  284. /* For convenience, set F'[k]=(k+F[k])/K */
  285. /* This saves some arithmetic in gsl_ran_discrete(); I find that
  286. * it doesn't actually make much difference.
  287. */
  288. for (k=0; k<Kevents; ++k) {
  289. (g->F)[k] += k;
  290. (g->F)[k] /= Kevents;
  291. }
  292. #endif
  293. free_stack(Bigs);
  294. free_stack(Smalls);
  295. free((char *)E);
  296. return g;
  297. }
  298. size_t
  299. gsl_ran_discrete(const gsl_rng *r, const gsl_ran_discrete_t *g)
  300. {
  301. size_t c=0;
  302. double u,f;
  303. u = gsl_rng_uniform(r);
  304. #if KNUTH_CONVENTION
  305. c = (u*(g->K));
  306. #else
  307. u *= g->K;
  308. c = u;
  309. u -= c;
  310. #endif
  311. f = (g->F)[c];
  312. /* fprintf(stderr,"c,f,u: %d %.4f %f\n",c,f,u); */
  313. if (f == 1.0) return c;
  314. if (u < f) {
  315. return c;
  316. }
  317. else {
  318. return (g->A)[c];
  319. }
  320. }
  321. void gsl_ran_discrete_free(gsl_ran_discrete_t *g)
  322. {
  323. free((char *)(g->A));
  324. free((char *)(g->F));
  325. free((char *)g);
  326. }
  327. double
  328. gsl_ran_discrete_pdf(size_t k, const gsl_ran_discrete_t *g)
  329. {
  330. size_t i,K;
  331. double f,p=0;
  332. K= g->K;
  333. if (k>K) return 0;
  334. for (i=0; i<K; ++i) {
  335. f = (g->F)[i];
  336. #if KNUTH_CONVENTION
  337. f = K*f-i;
  338. #endif
  339. if (i==k) {
  340. p += f;
  341. } else if (k == (g->A)[i]) {
  342. p += 1.0 - f;
  343. }
  344. }
  345. return p/K;
  346. }