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