Cleanup of old attributes. And addition of "top_tool" which is necessary for -clean...
[PsN.git] / modules / Math-Random / randlib.c
blobaac35b30c71d0afa6d2d77531e4dfb9744b82e22
1 #include "randlib.h"
2 #include <stdio.h>
3 #include <math.h>
4 #include <stdlib.h>
5 #define ABS(x) ((x) >= 0 ? (x) : -(x))
6 #define min(a,b) ((a) <= (b) ? (a) : (b))
7 #define max(a,b) ((a) >= (b) ? (a) : (b))
8 void ftnstop(char*);
10 double genbet(double aa,double bb)
12 **********************************************************************
13 double genbet(double aa,double bb)
14 GeNerate BETa random deviate
15 Function
16 Returns a single random deviate from the beta distribution with
17 parameters A and B. The density of the beta is
18 x^(a-1) * (1-x)^(b-1) / B(a,b) for 0 < x < 1
19 Arguments
20 aa --> First parameter of the beta distribution
22 bb --> Second parameter of the beta distribution
24 Method
25 R. C. H. Cheng
26 Generating Beta Variates with Nonintegral Shape Parameters
27 Communications of the ACM, 21:317-322 (1978)
28 (Algorithms BB and BC)
29 **********************************************************************
32 /* JJV changed expmax (log(1.0E38)==87.49823), and added minlog */
33 #define expmax 87.4982335337737
34 #define infnty 1.0E38
35 #define minlog 1.0E-37
36 static double olda = -1.0E37;
37 static double oldb = -1.0E37;
38 static double genbet,a,alpha,b,beta,delta,gamma,k1,k2,r,s,t,u1,u2,v,w,y,z;
39 static long qsame;
41 qsame = olda == aa && oldb == bb;
42 if(qsame) goto S20;
43 if(!(aa < minlog || bb < minlog)) goto S10;
44 fputs(" AA or BB < 1.0E-37 in GENBET - Abort!\n",stderr);
45 fprintf(stderr," AA: %16.6E BB %16.6E\n",aa,bb);
46 exit(1);
47 S10:
48 olda = aa;
49 oldb = bb;
50 S20:
51 if(!(min(aa,bb) > 1.0)) goto S100;
53 Algorithm BB
54 Initialize
56 if(qsame) goto S30;
57 a = min(aa,bb);
58 b = max(aa,bb);
59 alpha = a+b;
60 beta = sqrt((alpha-2.0)/(2.0*a*b-alpha));
61 gamma = a+1.0/beta;
62 S30:
63 u1 = ranf();
65 Step 1
67 u2 = ranf();
68 v = beta*log(u1/(1.0-u1));
69 /* JJV altered this */
70 if(v > expmax) goto S55;
72 * JJV added checker to see if a*exp(v) will overflow
73 * JJV S50 _was_ w = a*exp(v); also note here a > 1.0
75 w = exp(v);
76 if(w > infnty/a) goto S55;
77 w *= a;
78 goto S60;
79 S55:
80 w = infnty;
81 S60:
82 z = pow(u1,2.0)*u2;
83 r = gamma*v-1.38629436111989;
84 s = a+r-w;
86 Step 2
88 if(s+2.60943791243410 >= 5.0*z) goto S70;
90 Step 3
92 t = log(z);
93 if(s > t) goto S70;
95 * Step 4
97 * JJV added checker to see if log(alpha/(b+w)) will
98 * JJV overflow. If so, we count the log as -INF, and
99 * JJV consequently evaluate conditional as true, i.e.
100 * JJV the algorithm rejects the trial and starts over
101 * JJV May not need this here since alpha > 2.0
103 if(alpha/(b+w) < minlog) goto S30;
104 if(r+alpha*log(alpha/(b+w)) < t) goto S30;
105 S70:
107 Step 5
109 if(aa == a) {
110 genbet = w/(b+w);
111 } else {
112 genbet = b/(b+w);
114 goto S230;
115 S100:
117 Algorithm BC
118 Initialize
120 if(qsame) goto S110;
121 a = max(aa,bb);
122 b = min(aa,bb);
123 alpha = a+b;
124 beta = 1.0/b;
125 delta = 1.0+a-b;
126 k1 = delta*(1.38888888888889E-2+4.16666666666667E-2*b) /
127 (a*beta-0.777777777777778);
128 k2 = 0.25+(0.5+0.25/delta)*b;
129 S110:
130 S120:
131 u1 = ranf();
133 Step 1
135 u2 = ranf();
136 if(u1 >= 0.5) goto S130;
138 Step 2
140 y = u1*u2;
141 z = u1*y;
142 if(0.25*u2+z-y >= k1) goto S120;
143 goto S170;
144 S130:
146 Step 3
148 z = pow(u1,2.0)*u2;
149 if(!(z <= 0.25)) goto S160;
150 v = beta*log(u1/(1.0-u1));
152 * JJV instead of checking v > expmax at top, I will check
153 * JJV if a < 1, then check the appropriate values
155 if(a > 1.0) goto S135;
156 /* JJV a < 1 so it can help out if exp(v) would overflow */
157 if(v > expmax) goto S132;
158 w = a*exp(v);
159 goto S200;
160 S132:
161 w = v + log(a);
162 if(w > expmax) goto S140;
163 w = exp(w);
164 goto S200;
165 S135:
166 /* JJV in this case a > 1 */
167 if(v > expmax) goto S140;
168 w = exp(v);
169 if(w > infnty/a) goto S140;
170 w *= a;
171 goto S200;
172 S140:
173 w = infnty;
174 goto S200;
176 * JJV old code
177 * if(!(v > expmax)) goto S140;
178 * w = infnty;
179 * goto S150;
180 *S140:
181 * w = a*exp(v);
182 *S150:
183 * goto S200;
185 S160:
186 if(z >= k2) goto S120;
187 S170:
189 Step 4
190 Step 5
192 v = beta*log(u1/(1.0-u1));
193 /* JJV same kind of checking as above */
194 if(a > 1.0) goto S175;
195 /* JJV a < 1 so it can help out if exp(v) would overflow */
196 if(v > expmax) goto S172;
197 w = a*exp(v);
198 goto S190;
199 S172:
200 w = v + log(a);
201 if(w > expmax) goto S180;
202 w = exp(w);
203 goto S190;
204 S175:
205 /* JJV in this case a > 1.0 */
206 if(v > expmax) goto S180;
207 w = exp(v);
208 if(w > infnty/a) goto S180;
209 w *= a;
210 goto S190;
211 S180:
212 w = infnty;
214 * JJV old code
215 * if(!(v > expmax)) goto S180;
216 * w = infnty;
217 * goto S190;
218 *S180:
219 * w = a*exp(v);
221 S190:
223 * JJV here we also check to see if log overlows; if so, we treat it
224 * JJV as -INF, which means condition is true, i.e. restart
226 if(alpha/(b+w) < minlog) goto S120;
227 if(alpha*(log(alpha/(b+w))+v)-1.38629436111989 < log(z)) goto S120;
228 S200:
230 Step 6
232 if(a == aa) {
233 genbet = w/(b+w);
234 } else {
235 genbet = b/(b+w);
237 S230:
238 return genbet;
239 #undef expmax
240 #undef infnty
241 #undef minlog
244 double genchi(double df)
246 **********************************************************************
247 double genchi(double df)
248 Generate random value of CHIsquare variable
249 Function
250 Generates random deviate from the distribution of a chisquare
251 with DF degrees of freedom random variable.
252 Arguments
253 df --> Degrees of freedom of the chisquare
254 (Must be positive)
256 Method
257 Uses relation between chisquare and gamma.
258 **********************************************************************
261 static double genchi;
263 if(!(df <= 0.0)) goto S10;
264 fputs(" DF <= 0 in GENCHI - ABORT\n",stderr);
265 fprintf(stderr," Value of DF: %16.6E\n",df);
266 exit(1);
267 S10:
269 * JJV changed the code to call SGAMMA directly
270 * genchi = 2.0*gengam(1.0,df/2.0); <- OLD
272 genchi = 2.0*sgamma(df/2.0);
273 return genchi;
276 double genexp(double av)
278 **********************************************************************
279 double genexp(double av)
280 GENerate EXPonential random deviate
281 Function
282 Generates a single random deviate from an exponential
283 distribution with mean AV.
284 Arguments
285 av --> The mean of the exponential distribution from which
286 a random deviate is to be generated.
287 JJV (av >= 0)
288 Method
289 Renames SEXPO from TOMS as slightly modified by BWB to use RANF
290 instead of SUNIF.
291 For details see:
292 Ahrens, J.H. and Dieter, U.
293 Computer Methods for Sampling From the
294 Exponential and Normal Distributions.
295 Comm. ACM, 15,10 (Oct. 1972), 873 - 882.
296 **********************************************************************
299 static double genexp;
301 /* JJV added check that av >= 0 */
302 if(av >= 0.0) goto S10;
303 fputs(" AV < 0 in GENEXP - ABORT\n",stderr);
304 fprintf(stderr," Value of AV: %16.6E\n",av);
305 exit(1);
306 S10:
307 genexp = sexpo()*av;
308 return genexp;
311 double genf(double dfn,double dfd)
313 **********************************************************************
314 double genf(double dfn,double dfd)
315 GENerate random deviate from the F distribution
316 Function
317 Generates a random deviate from the F (variance ratio)
318 distribution with DFN degrees of freedom in the numerator
319 and DFD degrees of freedom in the denominator.
320 Arguments
321 dfn --> Numerator degrees of freedom
322 (Must be positive)
323 dfd --> Denominator degrees of freedom
324 (Must be positive)
325 Method
326 Directly generates ratio of chisquare variates
327 **********************************************************************
330 static double genf,xden,xnum;
332 if(!(dfn <= 0.0 || dfd <= 0.0)) goto S10;
333 fputs(" Degrees of freedom nonpositive in GENF - abort!\n",stderr);
334 fprintf(stderr," DFN value: %16.6E DFD value: %16.6E\n",dfn,dfd);
335 exit(1);
336 S10:
338 * JJV changed this to call SGAMMA directly
340 * GENF = ( GENCHI( DFN ) / DFN ) / ( GENCHI( DFD ) / DFD )
341 * xnum = genchi(dfn)/dfn; <- OLD
342 * xden = genchi(dfd)/dfd; <- OLD
344 xnum = 2.0*sgamma(dfn/2.0)/dfn;
345 xden = 2.0*sgamma(dfd/2.0)/dfd;
347 * JJV changed constant to prevent underflow at compile time.
348 * if(!(xden <= 9.999999999998E-39*xnum)) goto S20;
350 if(!(xden <= 1.0E-37*xnum)) goto S20;
351 fputs(" GENF - generated numbers would cause overflow\n",stderr);
352 fprintf(stderr," Numerator %16.6E Denominator %16.6E\n",xnum,xden);
354 * JJV changed next 2 lines to reflect constant change above in the
355 * JJV truncated value returned.
356 * fputs(" GENF returning 1.0E38\n",stderr);
357 * genf = 1.0E38;
359 fputs(" GENF returning 1.0E37\n",stderr);
360 genf = 1.0E37;
361 goto S30;
362 S20:
363 genf = xnum/xden;
364 S30:
365 return genf;
368 double gengam(double a,double r)
370 **********************************************************************
371 double gengam(double a,double r)
372 GENerates random deviates from GAMma distribution
373 Function
374 Generates random deviates from the gamma distribution whose
375 density is
376 (A**R)/Gamma(R) * X**(R-1) * Exp(-A*X)
377 Arguments
378 a --> Location parameter of Gamma distribution
379 JJV (a > 0)
380 r --> Shape parameter of Gamma distribution
381 JJV (r > 0)
382 Method
383 Renames SGAMMA from TOMS as slightly modified by BWB to use RANF
384 instead of SUNIF.
385 For details see:
386 (Case R >= 1.0)
387 Ahrens, J.H. and Dieter, U.
388 Generating Gamma Variates by a
389 Modified Rejection Technique.
390 Comm. ACM, 25,1 (Jan. 1982), 47 - 54.
391 Algorithm GD
392 JJV altered following to reflect argument ranges
393 (Case 0.0 < R < 1.0)
394 Ahrens, J.H. and Dieter, U.
395 Computer Methods for Sampling from Gamma,
396 Beta, Poisson and Binomial Distributions.
397 Computing, 12 (1974), 223-246/
398 Adapted algorithm GS.
399 **********************************************************************
402 static double gengam;
403 /* JJV added argument checker */
404 if(a > 0.0 && r > 0.0) goto S10;
405 fputs(" A or R nonpositive in GENGAM - abort!\n",stderr);
406 fprintf(stderr," A value: %16.6E R value: %16.6E\n",a,r);
407 exit(1);
408 S10:
409 gengam = sgamma(r);
410 gengam /= a;
411 return gengam;
414 void genmn(double *parm,double *x,double *work)
416 **********************************************************************
417 void genmn(double *parm,double *x,double *work)
418 GENerate Multivariate Normal random deviate
419 Arguments
420 parm --> Parameters needed to generate multivariate normal
421 deviates (MEANV and Cholesky decomposition of
422 COVM). Set by a previous call to SETGMN.
423 1 : 1 - size of deviate, P
424 2 : P + 1 - mean vector
425 P+2 : P*(P+3)/2 + 1 - upper half of cholesky
426 decomposition of cov matrix
427 x <-- Vector deviate generated.
428 work <--> Scratch array
429 Method
430 1) Generate P independent standard normal deviates - Ei ~ N(0,1)
431 2) Using Cholesky decomposition find A s.t. trans(A)*A = COVM
432 3) trans(A)E + MEANV ~ N(MEANV,COVM)
433 **********************************************************************
436 static long i,icount,j,p,D1,D2,D3,D4;
437 static double ae;
439 p = (long) (*parm);
441 Generate P independent normal deviates - WORK ~ N(0,1)
443 for(i=1; i<=p; i++) *(work+i-1) = snorm();
444 for(i=1,D3=1,D4=(p-i+D3)/D3; D4>0; D4--,i+=D3) {
446 PARM (P+2 : P*(P+3)/2 + 1) contains A, the Cholesky
447 decomposition of the desired covariance matrix.
448 trans(A)(1,1) = PARM(P+2)
449 trans(A)(2,1) = PARM(P+3)
450 trans(A)(2,2) = PARM(P+2+P)
451 trans(A)(3,1) = PARM(P+4)
452 trans(A)(3,2) = PARM(P+3+P)
453 trans(A)(3,3) = PARM(P+2-1+2P) ...
454 trans(A)*WORK + MEANV ~ N(MEANV,COVM)
456 icount = 0;
457 ae = 0.0;
458 for(j=1,D1=1,D2=(i-j+D1)/D1; D2>0; D2--,j+=D1) {
459 icount += (j-1);
460 ae += (*(parm+i+(j-1)*p-icount+p)**(work+j-1));
462 *(x+i-1) = ae+*(parm+i);
466 void genmul(long n,double *p,long ncat,long *ix)
468 **********************************************************************
470 void genmul(int n,double *p,int ncat,int *ix)
471 GENerate an observation from the MULtinomial distribution
472 Arguments
473 N --> Number of events that will be classified into one of
474 the categories 1..NCAT
475 P --> Vector of probabilities. P(i) is the probability that
476 an event will be classified into category i. Thus, P(i)
477 must be [0,1]. Only the first NCAT-1 P(i) must be defined
478 since P(NCAT) is 1.0 minus the sum of the first
479 NCAT-1 P(i).
480 NCAT --> Number of categories. Length of P and IX.
481 IX <-- Observation from multinomial distribution. All IX(i)
482 will be nonnegative and their sum will be N.
483 Method
484 Algorithm from page 559 of
486 Devroye, Luc
488 Non-Uniform Random Variate Generation. Springer-Verlag,
489 New York, 1986.
491 **********************************************************************
494 static double prob,ptot,sum;
495 static long i,icat,ntot;
496 if(n < 0) ftnstop("N < 0 in GENMUL");
497 if(ncat <= 1) ftnstop("NCAT <= 1 in GENMUL");
498 ptot = 0.0F;
499 for(i=0; i<ncat-1; i++) {
500 if(*(p+i) < 0.0F) ftnstop("Some P(i) < 0 in GENMUL");
501 if(*(p+i) > 1.0F) ftnstop("Some P(i) > 1 in GENMUL");
502 ptot += *(p+i);
504 if(ptot > 0.99999F) ftnstop("Sum of P(i) > 1 in GENMUL");
506 Initialize variables
508 ntot = n;
509 sum = 1.0F;
510 for(i=0; i<ncat; i++) ix[i] = 0;
512 Generate the observation
514 for(icat=0; icat<ncat-1; icat++) {
515 prob = *(p+icat)/sum;
516 *(ix+icat) = ignbin(ntot,prob);
517 ntot -= *(ix+icat);
518 if(ntot <= 0) return;
519 sum -= *(p+icat);
521 *(ix+ncat-1) = ntot;
523 Finished
525 return;
528 double gennch(double df,double xnonc)
530 **********************************************************************
531 double gennch(double df,double xnonc)
532 Generate random value of Noncentral CHIsquare variable
533 Function
534 Generates random deviate from the distribution of a noncentral
535 chisquare with DF degrees of freedom and noncentrality parameter
536 xnonc.
537 Arguments
538 df --> Degrees of freedom of the chisquare
539 (Must be >= 1.0)
540 xnonc --> Noncentrality parameter of the chisquare
541 (Must be >= 0.0)
542 Method
543 Uses fact that noncentral chisquare is the sum of a chisquare
544 deviate with DF-1 degrees of freedom plus the square of a normal
545 deviate with mean XNONC and standard deviation 1.
546 **********************************************************************
549 static double gennch;
551 if(!(df < 1.0 || xnonc < 0.0)) goto S10;
552 fputs("DF < 1 or XNONC < 0 in GENNCH - ABORT\n",stderr);
553 fprintf(stderr,"Value of DF: %16.6E Value of XNONC: %16.6E\n",df,xnonc);
554 exit(1);
555 /* JJV changed code to call SGAMMA, SNORM directly */
556 S10:
557 if(df >= 1.000000001) goto S20;
559 * JJV case df == 1.0
560 * gennch = pow(gennor(sqrt(xnonc),1.0),2.0); <- OLD
562 gennch = pow(snorm()+sqrt(xnonc),2.0);
563 goto S30;
564 S20:
566 * JJV case df > 1.0
567 * gennch = genchi(df-1.0)+pow(gennor(sqrt(xnonc),1.0),2.0); <- OLD
569 gennch = 2.0*sgamma((df-1.0)/2.0)+pow(snorm()+sqrt(xnonc),2.0);
570 S30:
571 return gennch;
574 double gennf(double dfn,double dfd,double xnonc)
576 **********************************************************************
577 double gennf(double dfn,double dfd,double xnonc)
578 GENerate random deviate from the Noncentral F distribution
579 Function
580 Generates a random deviate from the noncentral F (variance ratio)
581 distribution with DFN degrees of freedom in the numerator, and DFD
582 degrees of freedom in the denominator, and noncentrality parameter
583 XNONC.
584 Arguments
585 dfn --> Numerator degrees of freedom
586 (Must be >= 1.0)
587 dfd --> Denominator degrees of freedom
588 (Must be positive)
589 xnonc --> Noncentrality parameter
590 (Must be nonnegative)
591 Method
592 Directly generates ratio of noncentral numerator chisquare variate
593 to central denominator chisquare variate.
594 **********************************************************************
597 static double gennf,xden,xnum;
598 static long qcond;
600 /* JJV changed qcond, error message to allow dfn == 1.0 */
601 qcond = dfn < 1.0 || dfd <= 0.0 || xnonc < 0.0;
602 if(!qcond) goto S10;
603 fputs("In GENNF - Either (1) Numerator DF < 1.0 or\n",stderr);
604 fputs(" (2) Denominator DF <= 0.0 or\n",stderr);
605 fputs(" (3) Noncentrality parameter < 0.0\n",stderr);
606 fprintf(stderr,
607 "DFN value: %16.6E DFD value: %16.6E XNONC value: \n%16.6E\n",dfn,dfd,
608 xnonc);
609 exit(1);
610 S10:
612 * JJV changed the code to call SGAMMA and SNORM directly
613 * GENNF = ( GENNCH( DFN, XNONC ) / DFN ) / ( GENCHI( DFD ) / DFD )
614 * xnum = gennch(dfn,xnonc)/dfn; <- OLD
615 * xden = genchi(dfd)/dfd; <- OLD
617 if(dfn >= 1.000001) goto S20;
618 /* JJV case dfn == 1.0, dfn is counted as exactly 1.0 */
619 xnum = pow(snorm()+sqrt(xnonc),2.0);
620 goto S30;
621 S20:
622 /* JJV case df > 1.0 */
623 xnum = (2.0*sgamma((dfn-1.0)/2.0)+pow(snorm()+sqrt(xnonc),2.0))/dfn;
624 S30:
625 xden = 2.0*sgamma(dfd/2.0)/dfd;
627 * JJV changed constant to prevent underflow at compile time.
628 * if(!(xden <= 9.999999999998E-39*xnum)) goto S40;
630 if(!(xden <= 1.0E-37*xnum)) goto S40;
631 fputs(" GENNF - generated numbers would cause overflow\n",stderr);
632 fprintf(stderr," Numerator %16.6E Denominator %16.6E\n",xnum,xden);
634 * JJV changed next 2 lines to reflect constant change above in the
635 * JJV truncated value returned.
636 * fputs(" GENNF returning 1.0E38\n",stderr);
637 * gennf = 1.0E38;
639 fputs(" GENNF returning 1.0E37\n",stderr);
640 gennf = 1.0E37;
641 goto S50;
642 S40:
643 gennf = xnum/xden;
644 S50:
645 return gennf;
648 double gennor(double av,double sd)
650 **********************************************************************
651 double gennor(double av,double sd)
652 GENerate random deviate from a NORmal distribution
653 Function
654 Generates a single random deviate from a normal distribution
655 with mean, AV, and standard deviation, SD.
656 Arguments
657 av --> Mean of the normal distribution.
658 sd --> Standard deviation of the normal distribution.
659 JJV (sd >= 0)
660 Method
661 Renames SNORM from TOMS as slightly modified by BWB to use RANF
662 instead of SUNIF.
663 For details see:
664 Ahrens, J.H. and Dieter, U.
665 Extensions of Forsythe's Method for Random
666 Sampling from the Normal Distribution.
667 Math. Comput., 27,124 (Oct. 1973), 927 - 937.
668 **********************************************************************
671 static double gennor;
673 /* JJV added argument checker */
674 if(sd >= 0.0) goto S10;
675 fputs(" SD < 0 in GENNOR - ABORT\n",stderr);
676 fprintf(stderr," Value of SD: %16.6E\n",sd);
677 exit(1);
678 S10:
679 gennor = sd*snorm()+av;
680 return gennor;
683 void genprm(long *iarray,int larray)
685 **********************************************************************
686 void genprm(long *iarray,int larray)
687 GENerate random PeRMutation of iarray
688 Arguments
689 iarray <--> On output IARRAY is a random permutation of its
690 value on input
691 larray <--> Length of IARRAY
692 **********************************************************************
695 static long i,itmp,iwhich,D1,D2;
697 for(i=1,D1=1,D2=(larray-i+D1)/D1; D2>0; D2--,i+=D1) {
698 iwhich = ignuin(i,larray);
699 itmp = *(iarray+iwhich-1);
700 *(iarray+iwhich-1) = *(iarray+i-1);
701 *(iarray+i-1) = itmp;
705 double genunf(double low,double high)
707 **********************************************************************
708 double genunf(double low,double high)
709 GeNerate Uniform Real between LOW and HIGH
710 Function
711 Generates a real uniformly distributed between LOW and HIGH.
712 Arguments
713 low --> Low bound (exclusive) on real value to be generated
714 high --> High bound (exclusive) on real value to be generated
715 **********************************************************************
718 static double genunf;
720 if(!(low > high)) goto S10;
721 fprintf(stderr,"LOW > HIGH in GENUNF: LOW %16.6E HIGH: %16.6E\n",low,high);
722 fputs("Abort\n",stderr);
723 exit(1);
724 S10:
725 genunf = low+(high-low)*ranf();
726 return genunf;
729 void gscgn(long getset,long *g)
731 **********************************************************************
732 void gscgn(long getset,long *g)
733 Get/Set GeNerator
734 Gets or returns in G the number of the current generator
735 Arguments
736 getset --> 0 Get
737 1 Set
738 g <-- Number of the current random number generator (1..32)
739 **********************************************************************
742 #define numg 32L
743 static long curntg = 1;
744 if(getset == 0) *g = curntg;
745 else {
746 if(*g < 0 || *g > numg) {
747 fputs(" Generator number out of range in GSCGN\n",stderr);
748 exit(0);
750 curntg = *g;
752 #undef numg
755 void gsrgs(long getset,long *qvalue)
757 **********************************************************************
758 void gsrgs(long getset,long *qvalue)
759 Get/Set Random Generators Set
760 Gets or sets whether random generators set (initialized).
761 Initially (data statement) state is not set
762 If getset is 1 state is set to qvalue
763 If getset is 0 state returned in qvalue
764 **********************************************************************
767 static long qinit = 0;
769 if(getset == 0) *qvalue = qinit;
770 else qinit = *qvalue;
773 void gssst(long getset,long *qset)
775 **********************************************************************
776 void gssst(long getset,long *qset)
777 Get or Set whether Seed is Set
778 Initialize to Seed not Set
779 If getset is 1 sets state to Seed Set
780 If getset is 0 returns T in qset if Seed Set
781 Else returns F in qset
782 **********************************************************************
785 static long qstate = 0;
786 if(getset != 0) qstate = 1;
787 else *qset = qstate;
790 long ignbin(long n,double pp)
792 **********************************************************************
793 long ignbin(long n,double pp)
794 GENerate BINomial random deviate
795 Function
796 Generates a single random deviate from a binomial
797 distribution whose number of trials is N and whose
798 probability of an event in each trial is P.
799 Arguments
800 n --> The number of trials in the binomial distribution
801 from which a random deviate is to be generated.
802 JJV (N >= 0)
803 pp --> The probability of an event in each trial of the
804 binomial distribution from which a random deviate
805 is to be generated.
806 JJV (0.0 <= PP <= 1.0)
807 ignbin <-- A random deviate yielding the number of events
808 from N independent trials, each of which has
809 a probability of event P.
810 Method
811 This is algorithm BTPE from:
812 Kachitvichyanukul, V. and Schmeiser, B. W.
813 Binomial Random Variate Generation.
814 Communications of the ACM, 31, 2
815 (February, 1988) 216.
816 **********************************************************************
817 SUBROUTINE BTPEC(N,PP,ISEED,JX)
818 BINOMIAL RANDOM VARIATE GENERATOR
819 MEAN .LT. 30 -- INVERSE CDF
820 MEAN .GE. 30 -- ALGORITHM BTPE: ACCEPTANCE-REJECTION VIA
821 FOUR REGION COMPOSITION. THE FOUR REGIONS ARE A TRIANGLE
822 (SYMMETRIC IN THE CENTER), A PAIR OF PARALLELOGRAMS (ABOVE
823 THE TRIANGLE), AND EXPONENTIAL LEFT AND RIGHT TAILS.
824 BTPE REFERS TO BINOMIAL-TRIANGLE-PARALLELOGRAM-EXPONENTIAL.
825 BTPEC REFERS TO BTPE AND "COMBINED." THUS BTPE IS THE
826 RESEARCH AND BTPEC IS THE IMPLEMENTATION OF A COMPLETE
827 USABLE ALGORITHM.
828 REFERENCE: VORATAS KACHITVICHYANUKUL AND BRUCE SCHMEISER,
829 "BINOMIAL RANDOM VARIATE GENERATION,"
830 COMMUNICATIONS OF THE ACM, FORTHCOMING
831 WRITTEN: SEPTEMBER 1980.
832 LAST REVISED: MAY 1985, JULY 1987
833 REQUIRED SUBPROGRAM: RAND() -- A UNIFORM (0,1) RANDOM NUMBER
834 GENERATOR
835 ARGUMENTS
836 N : NUMBER OF BERNOULLI TRIALS (INPUT)
837 PP : PROBABILITY OF SUCCESS IN EACH TRIAL (INPUT)
838 ISEED: RANDOM NUMBER SEED (INPUT AND OUTPUT)
839 JX: RANDOMLY GENERATED OBSERVATION (OUTPUT)
840 VARIABLES
841 PSAVE: VALUE OF PP FROM THE LAST CALL TO BTPEC
842 NSAVE: VALUE OF N FROM THE LAST CALL TO BTPEC
843 XNP: VALUE OF THE MEAN FROM THE LAST CALL TO BTPEC
844 P: PROBABILITY USED IN THE GENERATION PHASE OF BTPEC
845 FFM: TEMPORARY VARIABLE EQUAL TO XNP + P
846 M: INTEGER VALUE OF THE CURRENT MODE
847 FM: FLOATING POINT VALUE OF THE CURRENT MODE
848 XNPQ: TEMPORARY VARIABLE USED IN SETUP AND SQUEEZING STEPS
849 P1: AREA OF THE TRIANGLE
850 C: HEIGHT OF THE PARALLELOGRAMS
851 XM: CENTER OF THE TRIANGLE
852 XL: LEFT END OF THE TRIANGLE
853 XR: RIGHT END OF THE TRIANGLE
854 AL: TEMPORARY VARIABLE
855 XLL: RATE FOR THE LEFT EXPONENTIAL TAIL
856 XLR: RATE FOR THE RIGHT EXPONENTIAL TAIL
857 P2: AREA OF THE PARALLELOGRAMS
858 P3: AREA OF THE LEFT EXPONENTIAL TAIL
859 P4: AREA OF THE RIGHT EXPONENTIAL TAIL
860 U: A U(0,P4) RANDOM VARIATE USED FIRST TO SELECT ONE OF THE
861 FOUR REGIONS AND THEN CONDITIONALLY TO GENERATE A VALUE
862 FROM THE REGION
863 V: A U(0,1) RANDOM NUMBER USED TO GENERATE THE RANDOM VALUE
864 (REGION 1) OR TRANSFORMED INTO THE VARIATE TO ACCEPT OR
865 REJECT THE CANDIDATE VALUE
866 IX: INTEGER CANDIDATE VALUE
867 X: PRELIMINARY CONTINUOUS CANDIDATE VALUE IN REGION 2 LOGIC
868 AND A FLOATING POINT IX IN THE ACCEPT/REJECT LOGIC
869 K: ABSOLUTE VALUE OF (IX-M)
870 F: THE HEIGHT OF THE SCALED DENSITY FUNCTION USED IN THE
871 ACCEPT/REJECT DECISION WHEN BOTH M AND IX ARE SMALL
872 ALSO USED IN THE INVERSE TRANSFORMATION
873 R: THE RATIO P/Q
874 G: CONSTANT USED IN CALCULATION OF PROBABILITY
875 MP: MODE PLUS ONE, THE LOWER INDEX FOR EXPLICIT CALCULATION
876 OF F WHEN IX IS GREATER THAN M
877 IX1: CANDIDATE VALUE PLUS ONE, THE LOWER INDEX FOR EXPLICIT
878 CALCULATION OF F WHEN IX IS LESS THAN M
879 I: INDEX FOR EXPLICIT CALCULATION OF F FOR BTPE
880 AMAXP: MAXIMUM ERROR OF THE LOGARITHM OF NORMAL BOUND
881 YNORM: LOGARITHM OF NORMAL BOUND
882 ALV: NATURAL LOGARITHM OF THE ACCEPT/REJECT VARIATE V
883 X1,F1,Z,W,Z2,X2,F2, AND W2 ARE TEMPORARY VARIABLES TO BE
884 USED IN THE FINAL ACCEPT/REJECT TEST
885 QN: PROBABILITY OF NO SUCCESS IN N TRIALS
886 REMARK
887 IX AND JX COULD LOGICALLY BE THE SAME VARIABLE, WHICH WOULD
888 SAVE A MEMORY POSITION AND A LINE OF CODE. HOWEVER, SOME
889 COMPILERS (E.G.,CDC MNF) OPTIMIZE BETTER WHEN THE ARGUMENTS
890 ARE NOT INVOLVED.
891 ISEED NEEDS TO BE DOUBLE PRECISION IF THE IMSL ROUTINE
892 GGUBFS IS USED TO GENERATE UNIFORM RANDOM NUMBER, OTHERWISE
893 TYPE OF ISEED SHOULD BE DICTATED BY THE UNIFORM GENERATOR
894 **********************************************************************
895 *****DETERMINE APPROPRIATE ALGORITHM AND WHETHER SETUP IS NECESSARY
898 /* JJV changed initial values to ridiculous values */
899 static double psave = -1.0E37;
900 static long nsave = -214748365;
901 static long ignbin,i,ix,ix1,k,m,mp,T1;
902 static double al,alv,amaxp,c,f,f1,f2,ffm,fm,g,p,p1,p2,p3,p4,q,qn,r,u,v,w,w2,x,x1,
903 x2,xl,xll,xlr,xm,xnp,xnpq,xr,ynorm,z,z2;
905 if(pp != psave) goto S10;
906 if(n != nsave) goto S20;
907 if(xnp < 30.0) goto S150;
908 goto S30;
909 S10:
911 *****SETUP, PERFORM ONLY WHEN PARAMETERS CHANGE
912 JJV added checks to ensure 0.0 <= PP <= 1.0
914 if(pp < 0.0F) ftnstop("PP < 0.0 in IGNBIN");
915 if(pp > 1.0F) ftnstop("PP > 1.0 in IGNBIN");
916 psave = pp;
917 p = min(psave,1.0-psave);
918 q = 1.0-p;
919 S20:
921 JJV added check to ensure N >= 0
923 if(n < 0L) ftnstop("N < 0 in IGNBIN");
924 xnp = n*p;
925 nsave = n;
926 if(xnp < 30.0) goto S140;
927 ffm = xnp+p;
928 m = ffm;
929 fm = m;
930 xnpq = xnp*q;
931 p1 = (long) (2.195*sqrt(xnpq)-4.6*q)+0.5;
932 xm = fm+0.5;
933 xl = xm-p1;
934 xr = xm+p1;
935 c = 0.134+20.5/(15.3+fm);
936 al = (ffm-xl)/(ffm-xl*p);
937 xll = al*(1.0+0.5*al);
938 al = (xr-ffm)/(xr*q);
939 xlr = al*(1.0+0.5*al);
940 p2 = p1*(1.0+c+c);
941 p3 = p2+c/xll;
942 p4 = p3+c/xlr;
943 S30:
945 *****GENERATE VARIATE
947 u = ranf()*p4;
948 v = ranf();
950 TRIANGULAR REGION
952 if(u > p1) goto S40;
953 ix = xm-p1*v+u;
954 goto S170;
955 S40:
957 PARALLELOGRAM REGION
959 if(u > p2) goto S50;
960 x = xl+(u-p1)/c;
961 v = v*c+1.0-ABS(xm-x)/p1;
962 if(v > 1.0 || v <= 0.0) goto S30;
963 ix = x;
964 goto S70;
965 S50:
967 LEFT TAIL
969 if(u > p3) goto S60;
970 ix = xl+log(v)/xll;
971 if(ix < 0) goto S30;
972 v *= ((u-p2)*xll);
973 goto S70;
974 S60:
976 RIGHT TAIL
978 ix = xr-log(v)/xlr;
979 if(ix > n) goto S30;
980 v *= ((u-p3)*xlr);
981 S70:
983 *****DETERMINE APPROPRIATE WAY TO PERFORM ACCEPT/REJECT TEST
985 k = ABS(ix-m);
986 if(k > 20 && k < xnpq/2-1) goto S130;
988 EXPLICIT EVALUATION
990 f = 1.0;
991 r = p/q;
992 g = (n+1)*r;
993 T1 = m-ix;
994 if(T1 < 0) goto S80;
995 else if(T1 == 0) goto S120;
996 else goto S100;
997 S80:
998 mp = m+1;
999 for(i=mp; i<=ix; i++) f *= (g/i-r);
1000 goto S120;
1001 S100:
1002 ix1 = ix+1;
1003 for(i=ix1; i<=m; i++) f /= (g/i-r);
1004 S120:
1005 if(v <= f) goto S170;
1006 goto S30;
1007 S130:
1009 SQUEEZING USING UPPER AND LOWER BOUNDS ON ALOG(F(X))
1011 amaxp = k/xnpq*((k*(k/3.0+0.625)+0.1666666666666)/xnpq+0.5);
1012 ynorm = -(k*k/(2.0*xnpq));
1013 alv = log(v);
1014 if(alv < ynorm-amaxp) goto S170;
1015 if(alv > ynorm+amaxp) goto S30;
1017 STIRLING'S FORMULA TO MACHINE ACCURACY FOR
1018 THE FINAL ACCEPTANCE/REJECTION TEST
1020 x1 = ix+1.0;
1021 f1 = fm+1.0;
1022 z = n+1.0-fm;
1023 w = n-ix+1.0;
1024 z2 = z*z;
1025 x2 = x1*x1;
1026 f2 = f1*f1;
1027 w2 = w*w;
1028 if(alv <= xm*log(f1/x1)+(n-m+0.5)*log(z/w)+(ix-m)*log(w*p/(x1*q))+(13860.0-
1029 (462.0-(132.0-(99.0-140.0/f2)/f2)/f2)/f2)/f1/166320.0+(13860.0-(462.0-
1030 (132.0-(99.0-140.0/z2)/z2)/z2)/z2)/z/166320.0+(13860.0-(462.0-(132.0-
1031 (99.0-140.0/x2)/x2)/x2)/x2)/x1/166320.0+(13860.0-(462.0-(132.0-(99.0
1032 -140.0/w2)/w2)/w2)/w2)/w/166320.0) goto S170;
1033 goto S30;
1034 S140:
1036 INVERSE CDF LOGIC FOR MEAN LESS THAN 30
1038 qn = pow(q,(double)n);
1039 r = p/q;
1040 g = r*(n+1);
1041 S150:
1042 ix = 0;
1043 f = qn;
1044 u = ranf();
1045 S160:
1046 if(u < f) goto S170;
1047 if(ix > 110) goto S150;
1048 u -= f;
1049 ix += 1;
1050 f *= (g/ix-r);
1051 goto S160;
1052 S170:
1053 if(psave > 0.5) ix = n-ix;
1054 ignbin = ix;
1055 return ignbin;
1058 long ignnbn(long n,double p)
1060 **********************************************************************
1062 long ignnbn(long n,double p)
1063 GENerate Negative BiNomial random deviate
1064 Function
1065 Generates a single random deviate from a negative binomial
1066 distribution.
1067 Arguments
1068 N --> The number of trials in the negative binomial distribution
1069 from which a random deviate is to be generated.
1070 JJV (N > 0)
1071 P --> The probability of an event.
1072 JJV (0.0 < P < 1.0)
1073 Method
1074 Algorithm from page 480 of
1076 Devroye, Luc
1078 Non-Uniform Random Variate Generation. Springer-Verlag,
1079 New York, 1986.
1080 **********************************************************************
1083 static long ignnbn;
1084 static double y,a,r;
1087 .. Executable Statements ..
1090 Check Arguments
1092 if(n <= 0L) ftnstop("N <= 0 in IGNNBN");
1093 if(p <= 0.0F) ftnstop("P <= 0.0 in IGNNBN");
1094 if(p >= 1.0F) ftnstop("P >= 1.0 in IGNNBN");
1096 Generate Y, a random gamma (n,(1-p)/p) variable
1097 JJV Note: the above parametrization is consistent with Devroye,
1098 JJV but gamma (p/(1-p),n) is the equivalent in our code
1100 r = (double)n;
1101 a = p/(1.0F-p);
1103 * JJV changed this to call SGAMMA directly
1104 * y = gengam(a,r); <- OLD
1106 y = sgamma(r)/a;
1108 Generate a random Poisson(y) variable
1110 ignnbn = ignpoi(y);
1111 return ignnbn;
1114 long ignpoi(double mu)
1116 **********************************************************************
1117 long ignpoi(double mu)
1118 GENerate POIsson random deviate
1119 Function
1120 Generates a single random deviate from a Poisson
1121 distribution with mean MU.
1122 Arguments
1123 mu --> The mean of the Poisson distribution from which
1124 a random deviate is to be generated.
1125 (mu >= 0.0)
1126 ignpoi <-- The random deviate.
1127 Method
1128 Renames KPOIS from TOMS as slightly modified by BWB to use RANF
1129 instead of SUNIF.
1130 For details see:
1131 Ahrens, J.H. and Dieter, U.
1132 Computer Generation of Poisson Deviates
1133 From Modified Normal Distributions.
1134 ACM Trans. Math. Software, 8, 2
1135 (June 1982),163-179
1136 **********************************************************************
1137 **********************************************************************
1140 P O I S S O N DISTRIBUTION
1143 **********************************************************************
1144 **********************************************************************
1146 FOR DETAILS SEE:
1148 AHRENS, J.H. AND DIETER, U.
1149 COMPUTER GENERATION OF POISSON DEVIATES
1150 FROM MODIFIED NORMAL DISTRIBUTIONS.
1151 ACM TRANS. MATH. SOFTWARE, 8,2 (JUNE 1982), 163 - 179.
1153 (SLIGHTLY MODIFIED VERSION OF THE PROGRAM IN THE ABOVE ARTICLE)
1155 **********************************************************************
1156 INTEGER FUNCTION IGNPOI(IR,MU)
1157 INPUT: IR=CURRENT STATE OF BASIC RANDOM NUMBER GENERATOR
1158 MU=MEAN MU OF THE POISSON DISTRIBUTION
1159 OUTPUT: IGNPOI=SAMPLE FROM THE POISSON-(MU)-DISTRIBUTION
1160 MUPREV=PREVIOUS MU, MUOLD=MU AT LAST EXECUTION OF STEP P OR B.
1161 TABLES: COEFFICIENTS A0-A7 FOR STEP F. FACTORIALS FACT
1162 COEFFICIENTS A(K) - FOR PX = FK*V*V*SUM(A(K)*V**K)-DEL
1163 SEPARATION OF CASES A AND B
1166 extern double fsign( double num, double sign );
1167 static double a0 = -0.5;
1168 static double a1 = 0.3333333343;
1169 static double a2 = -0.2499998565;
1170 static double a3 = 0.1999997049;
1171 static double a4 = -0.1666848753;
1172 static double a5 = 0.1428833286;
1173 static double a6 = -0.1241963125;
1174 static double a7 = 0.1101687109;
1175 static double a8 = -0.1142650302;
1176 static double a9 = 0.1055093006;
1177 /* JJV changed the initial values of MUPREV and MUOLD */
1178 static double muold = -1.0E37;
1179 static double muprev = -1.0E37;
1180 static double fact[10] = {
1181 1.0,1.0,2.0,6.0,24.0,120.0,720.0,5040.0,40320.0,362880.0
1183 /* JJV added ll to the list, for Case A */
1184 static long ignpoi,j,k,kflag,l,ll,m;
1185 static double b1,b2,c,c0,c1,c2,c3,d,del,difmuk,e,fk,fx,fy,g,omega,p,p0,px,py,q,s,
1186 t,u,v,x,xx,pp[35];
1188 if(mu == muprev) goto S10;
1189 if(mu < 10.0) goto S120;
1191 C A S E A. (RECALCULATION OF S,D,LL IF MU HAS CHANGED)
1192 JJV changed l in Case A to ll
1194 muprev = mu;
1195 s = sqrt(mu);
1196 d = 6.0*mu*mu;
1198 THE POISSON PROBABILITIES PK EXCEED THE DISCRETE NORMAL
1199 PROBABILITIES FK WHENEVER K >= M(MU). LL=IFIX(MU-1.1484)
1200 IS AN UPPER BOUND TO M(MU) FOR ALL MU >= 10 .
1202 ll = (long) (mu-1.1484);
1203 S10:
1205 STEP N. NORMAL SAMPLE - SNORM(IR) FOR STANDARD NORMAL DEVIATE
1207 g = mu+s*snorm();
1208 if(g < 0.0) goto S20;
1209 ignpoi = (long) (g);
1211 STEP I. IMMEDIATE ACCEPTANCE IF IGNPOI IS LARGE ENOUGH
1213 if(ignpoi >= ll) return ignpoi;
1215 STEP S. SQUEEZE ACCEPTANCE - SUNIF(IR) FOR (0,1)-SAMPLE U
1217 fk = (double)ignpoi;
1218 difmuk = mu-fk;
1219 u = ranf();
1220 if(d*u >= difmuk*difmuk*difmuk) return ignpoi;
1221 S20:
1223 STEP P. PREPARATIONS FOR STEPS Q AND H.
1224 (RECALCULATIONS OF PARAMETERS IF NECESSARY)
1225 .3989423=(2*PI)**(-.5) .416667E-1=1./24. .1428571=1./7.
1226 THE QUANTITIES B1, B2, C3, C2, C1, C0 ARE FOR THE HERMITE
1227 APPROXIMATIONS TO THE DISCRETE NORMAL PROBABILITIES FK.
1228 C=.1069/MU GUARANTEES MAJORIZATION BY THE 'HAT'-FUNCTION.
1230 if(mu == muold) goto S30;
1231 muold = mu;
1232 omega = 0.398942280401433/s;
1233 b1 = 4.16666666666667E-2/mu;
1234 b2 = 0.3*b1*b1;
1235 c3 = 0.142857142857143*b1*b2;
1236 c2 = b2-15.0*c3;
1237 c1 = b1-6.0*b2+45.0*c3;
1238 c0 = 1.0-b1+3.0*b2-15.0*c3;
1239 c = 0.1069/mu;
1240 S30:
1241 if(g < 0.0) goto S50;
1243 'SUBROUTINE' F IS CALLED (KFLAG=0 FOR CORRECT RETURN)
1245 kflag = 0;
1246 goto S70;
1247 S40:
1249 STEP Q. QUOTIENT ACCEPTANCE (RARE CASE)
1251 if(fy-u*fy <= py*exp(px-fx)) return ignpoi;
1252 S50:
1254 STEP E. EXPONENTIAL SAMPLE - SEXPO(IR) FOR STANDARD EXPONENTIAL
1255 DEVIATE E AND SAMPLE T FROM THE LAPLACE 'HAT'
1256 (IF T <= -.6744 THEN PK < FK FOR ALL MU >= 10.)
1258 e = sexpo();
1259 u = ranf();
1260 u += (u-1.0);
1261 t = 1.8+fsign(e,u);
1262 if(t <= -0.6744) goto S50;
1263 ignpoi = (long) (mu+s*t);
1264 fk = (double)ignpoi;
1265 difmuk = mu-fk;
1267 'SUBROUTINE' F IS CALLED (KFLAG=1 FOR CORRECT RETURN)
1269 kflag = 1;
1270 goto S70;
1271 S60:
1273 STEP H. HAT ACCEPTANCE (E IS REPEATED ON REJECTION)
1275 if(c*fabs(u) > py*exp(px+e)-fy*exp(fx+e)) goto S50;
1276 return ignpoi;
1277 S70:
1279 STEP F. 'SUBROUTINE' F. CALCULATION OF PX,PY,FX,FY.
1280 CASE IGNPOI .LT. 10 USES FACTORIALS FROM TABLE FACT
1282 if(ignpoi >= 10) goto S80;
1283 px = -mu;
1284 py = pow(mu,(double)ignpoi)/ *(fact+ignpoi);
1285 goto S110;
1286 S80:
1288 CASE IGNPOI .GE. 10 USES POLYNOMIAL APPROXIMATION
1289 A0-A7 FOR ACCURACY WHEN ADVISABLE
1290 .8333333E-1=1./12. .3989423=(2*PI)**(-.5)
1292 del = 8.33333333E-2/fk;
1293 del -= (4.8*del*del*del);
1294 v = difmuk/fk;
1295 if(fabs(v) <= 0.25) goto S90;
1296 px = fk*log(1.0+v)-difmuk-del;
1297 goto S100;
1298 S90:
1299 px = fk*v*v*((((((((a8*v+a7)*v+a6)*v+a5)*v+a4)*v+a3)*v+a2)*v+a1)*v+a0)-del;
1300 S100:
1301 py = 0.398942280401433/sqrt(fk);
1302 S110:
1303 x = (0.5-difmuk)/s;
1304 xx = x*x;
1305 fx = -0.5*xx;
1306 fy = omega*(((c3*xx+c2)*xx+c1)*xx+c0);
1307 if(kflag <= 0) goto S40;
1308 goto S60;
1309 S120:
1311 C A S E B. (START NEW TABLE AND CALCULATE P0 IF NECESSARY)
1312 JJV changed MUPREV assignment to initial value
1314 muprev = -1.0E37;
1315 if(mu == muold) goto S130;
1316 /* JJV added argument checker here */
1317 if(mu >= 0.0) goto S125;
1318 fprintf(stderr,"MU < 0 in IGNPOI: MU %16.6E\n",mu);
1319 fputs("Abort\n",stderr);
1320 exit(1);
1321 S125:
1322 muold = mu;
1323 m = max(1L,(long) (mu));
1324 l = 0;
1325 p = exp(-mu);
1326 q = p0 = p;
1327 S130:
1329 STEP U. UNIFORM SAMPLE FOR INVERSION METHOD
1331 u = ranf();
1332 ignpoi = 0;
1333 if(u <= p0) return ignpoi;
1335 STEP T. TABLE COMPARISON UNTIL THE END PP(L) OF THE
1336 PP-TABLE OF CUMULATIVE POISSON PROBABILITIES
1337 (0.458=PP(9) FOR MU=10)
1339 if(l == 0) goto S150;
1340 j = 1;
1341 if(u > 0.458) j = min(l,m);
1342 for(k=j; k<=l; k++) {
1343 if(u <= *(pp+k-1)) goto S180;
1345 if(l == 35) goto S130;
1346 S150:
1348 STEP C. CREATION OF NEW POISSON PROBABILITIES P
1349 AND THEIR CUMULATIVES Q=PP(K)
1351 l += 1;
1352 for(k=l; k<=35; k++) {
1353 p = p*mu/(double)k;
1354 q += p;
1355 *(pp+k-1) = q;
1356 if(u <= q) goto S170;
1358 l = 35;
1359 goto S130;
1360 S170:
1361 l = k;
1362 S180:
1363 ignpoi = k;
1364 return ignpoi;
1367 long ignuin(long low,long high)
1369 **********************************************************************
1370 long ignuin(long low,long high)
1371 GeNerate Uniform INteger
1372 Function
1373 Generates an integer uniformly distributed between LOW and HIGH.
1374 Arguments
1375 low --> Low bound (inclusive) on integer value to be generated
1376 high --> High bound (inclusive) on integer value to be generated
1377 Note
1378 If (HIGH-LOW) > 2,147,483,561 prints error message on * unit and
1379 stops the program.
1380 **********************************************************************
1381 IGNLGI generates integers between 1 and 2147483562
1382 MAXNUM is 1 less than maximum generable value
1385 #define maxnum 2147483561L
1386 static long ignuin,ign,maxnow,range,ranp1;
1388 if(!(low > high)) goto S10;
1389 fputs(" low > high in ignuin - ABORT\n",stderr);
1390 exit(1);
1392 S10:
1393 range = high-low;
1394 if(!(range > maxnum)) goto S20;
1395 fputs(" high - low too large in ignuin - ABORT\n",stderr);
1396 exit(1);
1398 S20:
1399 if(!(low == high)) goto S30;
1400 ignuin = low;
1401 return ignuin;
1403 S30:
1405 Number to be generated should be in range 0..RANGE
1406 Set MAXNOW so that the number of integers in 0..MAXNOW is an
1407 integral multiple of the number in 0..RANGE
1409 ranp1 = range+1;
1410 maxnow = maxnum/ranp1*ranp1;
1411 S40:
1412 ign = ignlgi()-1;
1413 if(!(ign <= maxnow)) goto S40;
1414 ignuin = low+ign%ranp1;
1415 return ignuin;
1416 #undef maxnum
1417 #undef err1
1418 #undef err2
1421 long lennob( char *str )
1423 Returns the length of str ignoring trailing blanks but not
1424 other white space.
1427 long i, i_nb;
1429 for (i=0, i_nb= -1L; *(str+i); i++)
1430 if ( *(str+i) != ' ' ) i_nb = i;
1431 return (i_nb+1);
1434 long mltmod(long a,long s,long m)
1436 **********************************************************************
1437 long mltmod(long a,long s,long m)
1438 Returns (A*S) MOD M
1439 This is a transcription from Pascal to C of routine
1440 MultMod_Decompos from the paper
1441 L'Ecuyer, P. and Cote, S. "Implementing a Random Number Package
1442 with Splitting Facilities." ACM Transactions on Mathematical
1443 Software, 17:98-111 (1991)
1444 Arguments
1445 a, s, m -->
1446 WGR, 12/19/00: replaced S10, S20, etc. with C blocks {} per
1447 original paper.
1448 **********************************************************************
1451 #define h 32768L
1452 static long a0,a1,k,p,q,qh,rh;
1454 H = 2**((b-2)/2) where b = 32 because we are using a 32 bit
1455 machine. On a different machine recompute H.
1457 if (a <= 0 || a >= m || s <= 0 || s >= m) {
1458 fputs(" a, m, s out of order in mltmod - ABORT!\n",stderr);
1459 fprintf(stderr," a = %12ld s = %12ld m = %12ld\n",a,s,m);
1460 fputs(" mltmod requires: 0 < a < m; 0 < s < m\n",stderr);
1461 exit(1);
1464 if (a < h) {
1465 a0 = a;
1466 p = 0;
1467 } else {
1468 a1 = a/h;
1469 a0 = a - h*a1;
1470 qh = m/h;
1471 rh = m - h*qh;
1472 if (a1 >= h) { /* A2=1 */
1473 a1 -= h;
1474 k = s/qh;
1475 p = h*(s-k*qh) - k*rh;
1476 while (p < 0) { p += m; }
1477 } else {
1478 p = 0;
1481 P = (A2*S*H)MOD M
1483 if (a1 != 0) {
1484 q = m/a1;
1485 k = s/q;
1486 p -= k*(m - a1*q);
1487 if (p > 0) { p -= m; }
1488 p += a1*(s - k*q);
1489 while (p < 0) { p += m; }
1492 P = ((A2*H + A1)*S)MOD M
1494 k = p/qh;
1495 p = h*(p-k*qh) - k*rh;
1496 while (p < 0) { p += m; }
1499 P = ((A2*H + A1)*H*S)MOD M
1501 if (a0 != 0) {
1502 q = m/a0;
1503 k = s/q;
1504 p -= k*(m-a0*q);
1505 if (p > 0) { p -= m; }
1506 p += a0*(s-k*q);
1507 while (p < 0) { p += m; }
1509 return p;
1510 #undef h
1513 void phrtsd(char* phrase,long *seed1,long *seed2)
1515 **********************************************************************
1516 void phrtsd(char* phrase,long *seed1,long *seed2)
1517 PHRase To SeeDs
1519 Function
1521 Uses a phrase (character string) to generate two seeds for the RGN
1522 random number generator.
1523 Arguments
1524 phrase --> Phrase to be used for random number generation
1526 seed1 <-- First seed for generator
1528 seed2 <-- Second seed for generator
1530 Note
1532 Trailing blanks are eliminated before the seeds are generated.
1533 Generated seed values will fall in the range 1..2^30
1534 (1..1,073,741,824)
1535 **********************************************************************
1539 static char table[] =
1540 "abcdefghijklmnopqrstuvwxyz\
1541 ABCDEFGHIJKLMNOPQRSTUVWXYZ\
1542 0123456789\
1543 !@#$%^&*()_+[];:'\\\"<>?,./ "; /* WGR added space, 5/19/1999 */
1545 long ix;
1547 static long twop30 = 1073741824L;
1548 static long shift[5] = {
1549 1L,64L,4096L,262144L,16777216L
1551 static long i,ichr,j,lphr,values[5];
1552 extern long lennob(char *str);
1554 *seed1 = 1234567890L;
1555 *seed2 = 123456789L;
1556 lphr = lennob(phrase);
1557 if(lphr < 1) return;
1558 for(i=0; i<=(lphr-1); i++) {
1559 for (ix=0; table[ix]; ix++) if (*(phrase+i) == table[ix]) break;
1560 /* JJV added ix++; to bring index in line with fortran's index*/
1561 ix++;
1562 if (!table[ix]) ix = 0;
1563 ichr = ix % 64;
1564 if(ichr == 0) ichr = 63;
1565 for(j=1; j<=5; j++) {
1566 *(values+j-1) = ichr-j;
1567 if(*(values+j-1) < 1) *(values+j-1) += 63;
1569 for(j=1; j<=5; j++) {
1570 *seed1 = ( *seed1+*(shift+j-1)**(values+j-1) ) % twop30;
1571 *seed2 = ( *seed2+*(shift+j-1)**(values+6-j-1) ) % twop30;
1576 double ranf(void)
1578 **********************************************************************
1579 double ranf(void)
1580 RANDom number generator as a Function
1581 Returns a random floating point number from a uniform distribution
1582 over 0 - 1 (endpoints of this interval are not returned) using the
1583 current generator.
1584 This is a transcription from Pascal to C of routine
1585 Uniform_01 from the paper
1586 L'Ecuyer, P. and Cote, S. "Implementing a Random Number Package
1587 with Splitting Facilities." ACM Transactions on Mathematical
1588 Software, 17:98-111 (1991)
1589 WGR, 2/12/01: increased precision.
1590 **********************************************************************
1593 static double ranf;
1595 4.656613057E-10 is 1/M1 M1 is set in a data statement in IGNLGI
1596 and is currently 2147483563. If M1 changes, change this also.
1598 ranf = ignlgi()*4.65661305739177E-10;
1599 return ranf;
1602 void setgmn(double *meanv,double *covm,long p,double *parm)
1604 **********************************************************************
1605 void setgmn(double *meanv,double *covm,long p,double *parm)
1606 SET Generate Multivariate Normal random deviate
1607 Function
1608 Places P, MEANV, and the Cholesky factorization of COVM
1609 in GENMN.
1610 Arguments
1611 meanv --> Mean vector of multivariate normal distribution.
1612 covm <--> (Input) Covariance matrix of the multivariate
1613 normal distribution
1614 (Output) Destroyed on output
1615 p --> Dimension of the normal, or length of MEANV.
1616 parm <-- Array of parameters needed to generate multivariate norma
1617 deviates (P, MEANV and Cholesky decomposition of
1618 COVM).
1619 1 : 1 - P
1620 2 : P + 1 - MEANV
1621 P+2 : P*(P+3)/2 + 1 - Cholesky decomposition of COVM
1622 Needed dimension is (p*(p+3)/2 + 1)
1623 **********************************************************************
1626 extern void spofa(double *a,long lda,long n,long *info);
1627 static long T1;
1628 static long i,icount,info,j,D2,D3,D4,D5;
1629 T1 = p*(p+3)/2+1;
1631 TEST THE INPUT
1633 if(!(p <= 0)) goto S10;
1634 fputs("P nonpositive in SETGMN\n",stderr);
1635 fprintf(stderr,"Value of P: %12ld\n",p);
1636 exit(1);
1637 S10:
1638 *parm = p;
1640 PUT P AND MEANV INTO PARM
1642 for(i=2,D2=1,D3=(p+1-i+D2)/D2; D3>0; D3--,i+=D2) *(parm+i-1) = *(meanv+i-2);
1644 Cholesky decomposition to find A s.t. trans(A)*(A) = COVM
1646 spofa(covm,p,p,&info);
1647 if(!(info != 0)) goto S30;
1648 fputs(" COVM not positive definite in SETGMN\n",stderr);
1649 exit(1);
1650 S30:
1651 icount = p+1;
1653 PUT UPPER HALF OF A, WHICH IS NOW THE CHOLESKY FACTOR, INTO PARM
1654 COVM(1,1) = PARM(P+2)
1655 COVM(1,2) = PARM(P+3)
1657 COVM(1,P) = PARM(2P+1)
1658 COVM(2,2) = PARM(2P+2) ...
1660 for(i=1,D4=1,D5=(p-i+D4)/D4; D5>0; D5--,i+=D4) {
1661 for(j=i-1; j<p; j++) {
1662 icount += 1;
1663 *(parm+icount-1) = *(covm+i-1+j*p);
1668 double sexpo(void)
1670 **********************************************************************
1673 (STANDARD-) E X P O N E N T I A L DISTRIBUTION
1676 **********************************************************************
1677 **********************************************************************
1679 FOR DETAILS SEE:
1681 AHRENS, J.H. AND DIETER, U.
1682 COMPUTER METHODS FOR SAMPLING FROM THE
1683 EXPONENTIAL AND NORMAL DISTRIBUTIONS.
1684 COMM. ACM, 15,10 (OCT. 1972), 873 - 882.
1686 ALL STATEMENT NUMBERS CORRESPOND TO THE STEPS OF ALGORITHM
1687 'SA' IN THE ABOVE PAPER (SLIGHTLY MODIFIED IMPLEMENTATION)
1689 Modified by Barry W. Brown, Feb 3, 1988 to use RANF instead of
1690 SUNIF. The argument IR thus goes away.
1692 **********************************************************************
1693 Q(N) = SUM(ALOG(2.0)**K/K!) K=1,..,N , THE HIGHEST N
1694 (HERE 8) IS DETERMINED BY Q(N)=1.0 WITHIN STANDARD PRECISION
1697 static double q[8] = {
1698 0.69314718055995, 0.93337368751905, 0.98887779618387, 0.99849592529150,
1699 0.99982928110614, 0.99998331641007, 0.99999856914388, 0.99999989069256
1701 static long i;
1702 static double sexpo,a,u,ustar,umin;
1703 static double *q1 = q;
1704 a = 0.0;
1705 u = ranf();
1706 goto S30;
1707 S20:
1708 a += *q1;
1709 S30:
1710 u += u;
1712 * JJV changed the following to reflect the true algorithm and prevent
1713 * JJV unpredictable behavior if U is initially 0.5.
1714 * if(u <= 1.0) goto S20;
1716 if(u < 1.0) goto S20;
1717 u -= 1.0;
1718 if(u > *q1) goto S60;
1719 sexpo = a+u;
1720 return sexpo;
1721 S60:
1722 i = 1;
1723 ustar = ranf();
1724 umin = ustar;
1725 S70:
1726 ustar = ranf();
1727 if(ustar < umin) umin = ustar;
1728 i += 1;
1729 if(u > *(q+i-1)) goto S70;
1730 sexpo = a+umin**q1;
1731 return sexpo;
1734 double sgamma(double a)
1736 **********************************************************************
1739 (STANDARD-) G A M M A DISTRIBUTION
1742 **********************************************************************
1743 **********************************************************************
1745 PARAMETER A >= 1.0 !
1747 **********************************************************************
1749 FOR DETAILS SEE:
1751 AHRENS, J.H. AND DIETER, U.
1752 GENERATING GAMMA VARIATES BY A
1753 MODIFIED REJECTION TECHNIQUE.
1754 COMM. ACM, 25,1 (JAN. 1982), 47 - 54.
1756 STEP NUMBERS CORRESPOND TO ALGORITHM 'GD' IN THE ABOVE PAPER
1757 (STRAIGHTFORWARD IMPLEMENTATION)
1759 Modified by Barry W. Brown, Feb 3, 1988 to use RANF instead of
1760 SUNIF. The argument IR thus goes away.
1762 **********************************************************************
1764 PARAMETER 0.0 < A < 1.0 !
1766 **********************************************************************
1768 FOR DETAILS SEE:
1770 AHRENS, J.H. AND DIETER, U.
1771 COMPUTER METHODS FOR SAMPLING FROM GAMMA,
1772 BETA, POISSON AND BINOMIAL DISTRIBUTIONS.
1773 COMPUTING, 12 (1974), 223 - 246.
1775 (ADAPTED IMPLEMENTATION OF ALGORITHM 'GS' IN THE ABOVE PAPER)
1777 **********************************************************************
1778 INPUT: A =PARAMETER (MEAN) OF THE STANDARD GAMMA DISTRIBUTION
1779 OUTPUT: SGAMMA = SAMPLE FROM THE GAMMA-(A)-DISTRIBUTION
1780 COEFFICIENTS Q(K) - FOR Q0 = SUM(Q(K)*A**(-K))
1781 COEFFICIENTS A(K) - FOR Q = Q0+(T*T/2)*SUM(A(K)*V**K)
1782 COEFFICIENTS E(K) - FOR EXP(Q)-1 = SUM(E(K)*Q**K)
1783 PREVIOUS A PRE-SET TO ZERO - AA IS A', AAA IS A"
1784 SQRT32 IS THE SQUAREROOT OF 32 = 5.656854249492380
1787 extern double fsign( double num, double sign );
1788 static double q1 = 4.16666664E-2;
1789 static double q2 = 2.08333723E-2;
1790 static double q3 = 7.9849875E-3;
1791 static double q4 = 1.5746717E-3;
1792 static double q5 = -3.349403E-4;
1793 static double q6 = 3.340332E-4;
1794 static double q7 = 6.053049E-4;
1795 static double q8 = -4.701849E-4;
1796 static double q9 = 1.710320E-4;
1797 static double a1 = 0.333333333;
1798 static double a2 = -0.249999949;
1799 static double a3 = 0.199999867;
1800 static double a4 = -0.166677482;
1801 static double a5 = 0.142873973;
1802 static double a6 = -0.124385581;
1803 static double a7 = 0.110368310;
1804 static double a8 = -0.112750886;
1805 static double a9 = 0.104089866;
1806 static double e1 = 1.0;
1807 static double e2 = 0.499999994;
1808 static double e3 = 0.166666848;
1809 static double e4 = 4.1664508E-2;
1810 static double e5 = 8.345522E-3;
1811 static double e6 = 1.353826E-3;
1812 static double e7 = 2.47453E-4;
1813 static double aa = 0.0;
1814 static double aaa = 0.0;
1815 static double sqrt32 = 5.65685424949238;
1816 /* JJV added b0 to fix rare and subtle bug */
1817 static double sgamma,s2,s,d,t,x,u,r,q0,b,b0,si,c,v,q,e,w,p;
1818 if(a == aa) goto S10;
1819 if(a < 1.0) goto S120;
1821 STEP 1: RECALCULATIONS OF S2,S,D IF A HAS CHANGED
1823 aa = a;
1824 s2 = a-0.5;
1825 s = sqrt(s2);
1826 d = sqrt32-12.0*s;
1827 S10:
1829 STEP 2: T=STANDARD NORMAL DEVIATE,
1830 X=(S,1/2)-NORMAL DEVIATE.
1831 IMMEDIATE ACCEPTANCE (I)
1833 t = snorm();
1834 x = s+0.5*t;
1835 sgamma = x*x;
1836 if(t >= 0.0) return sgamma;
1838 STEP 3: U= 0,1 -UNIFORM SAMPLE. SQUEEZE ACCEPTANCE (S)
1840 u = ranf();
1841 if(d*u <= t*t*t) return sgamma;
1843 STEP 4: RECALCULATIONS OF Q0,B,SI,C IF NECESSARY
1845 if(a == aaa) goto S40;
1846 aaa = a;
1847 r = 1.0/a;
1848 q0 = ((((((((q9*r+q8)*r+q7)*r+q6)*r+q5)*r+q4)*r+q3)*r+q2)*r+q1)*r;
1850 APPROXIMATION DEPENDING ON SIZE OF PARAMETER A
1851 THE CONSTANTS IN THE EXPRESSIONS FOR B, SI AND
1852 C WERE ESTABLISHED BY NUMERICAL EXPERIMENTS
1854 if(a <= 3.686) goto S30;
1855 if(a <= 13.022) goto S20;
1857 CASE 3: A .GT. 13.022
1859 b = 1.77;
1860 si = 0.75;
1861 c = 0.1515/s;
1862 goto S40;
1863 S20:
1865 CASE 2: 3.686 .LT. A .LE. 13.022
1867 b = 1.654+7.6E-3*s2;
1868 si = 1.68/s+0.275;
1869 c = 6.2E-2/s+2.4E-2;
1870 goto S40;
1871 S30:
1873 CASE 1: A .LE. 3.686
1875 b = 0.463+s+0.178*s2;
1876 si = 1.235;
1877 c = 0.195/s-7.9E-2+1.6E-1*s;
1878 S40:
1880 STEP 5: NO QUOTIENT TEST IF X NOT POSITIVE
1882 if(x <= 0.0) goto S70;
1884 STEP 6: CALCULATION OF V AND QUOTIENT Q
1886 v = t/(s+s);
1887 if(fabs(v) <= 0.25) goto S50;
1888 q = q0-s*t+0.25*t*t+(s2+s2)*log(1.0+v);
1889 goto S60;
1890 S50:
1891 q = q0+0.5*t*t*((((((((a9*v+a8)*v+a7)*v+a6)*v+a5)*v+a4)*v+a3)*v+a2)*v+a1)*v;
1892 S60:
1894 STEP 7: QUOTIENT ACCEPTANCE (Q)
1896 if(log(1.0-u) <= q) return sgamma;
1897 S70:
1899 STEP 8: E=STANDARD EXPONENTIAL DEVIATE
1900 U= 0,1 -UNIFORM DEVIATE
1901 T=(B,SI)-DOUBLE EXPONENTIAL (LAPLACE) SAMPLE
1903 e = sexpo();
1904 u = ranf();
1905 u += (u-1.0);
1906 t = b+fsign(si*e,u);
1908 STEP 9: REJECTION IF T .LT. TAU(1) = -.71874483771719
1910 if(t < -0.71874483771719) goto S70;
1912 STEP 10: CALCULATION OF V AND QUOTIENT Q
1914 v = t/(s+s);
1915 if(fabs(v) <= 0.25) goto S80;
1916 q = q0-s*t+0.25*t*t+(s2+s2)*log(1.0+v);
1917 goto S90;
1918 S80:
1919 q = q0+0.5*t*t*((((((((a9*v+a8)*v+a7)*v+a6)*v+a5)*v+a4)*v+a3)*v+a2)*v+a1)*v;
1920 S90:
1922 STEP 11: HAT ACCEPTANCE (H) (IF Q NOT POSITIVE GO TO STEP 8)
1924 if(q <= 0.0) goto S70;
1925 if(q <= 0.5) goto S100;
1927 * JJV modified the code through line 115 to handle large Q case
1929 if(q < 15.0) goto S95;
1931 * JJV Here Q is large enough that Q = log(exp(Q) - 1.0) (for real Q)
1932 * JJV so reformulate test at 110 in terms of one EXP, if not too big
1933 * JJV 87.49823 is close to the largest real which can be
1934 * JJV exponentiated (87.49823 = log(1.0E38))
1936 if((q+e-0.5*t*t) > 87.4982335337737) goto S115;
1937 if(c*fabs(u) > exp(q+e-0.5*t*t)) goto S70;
1938 goto S115;
1939 S95:
1940 w = exp(q)-1.0;
1941 goto S110;
1942 S100:
1943 w = ((((((e7*q+e6)*q+e5)*q+e4)*q+e3)*q+e2)*q+e1)*q;
1944 S110:
1946 IF T IS REJECTED, SAMPLE AGAIN AT STEP 8
1948 if(c*fabs(u) > w*exp(e-0.5*t*t)) goto S70;
1949 S115:
1950 x = s+0.5*t;
1951 sgamma = x*x;
1952 return sgamma;
1953 S120:
1955 ALTERNATE METHOD FOR PARAMETERS A BELOW 1 (.3678794=EXP(-1.))
1957 JJV changed B to B0 (which was added to declarations for this)
1958 JJV in 120 to END to fix rare and subtle bug.
1959 JJV Line: 'aa = 0.0' was removed (unnecessary, wasteful).
1960 JJV Reasons: the state of AA only serves to tell the A >= 1.0
1961 JJV case if certain A-dependent constants need to be recalculated.
1962 JJV The A < 1.0 case (here) no longer changes any of these, and
1963 JJV the recalculation of B (which used to change with an
1964 JJV A < 1.0 call) is governed by the state of AAA anyway.
1965 aa = 0.0;
1967 b0 = 1.0+ 0.3678794411714423*a;
1968 S130:
1969 p = b0*ranf();
1970 if(p >= 1.0) goto S140;
1971 sgamma = exp(log(p)/ a);
1972 if(sexpo() < sgamma) goto S130;
1973 return sgamma;
1974 S140:
1975 sgamma = -log((b0-p)/ a);
1976 if(sexpo() < (1.0-a)*log(sgamma)) goto S130;
1977 return sgamma;
1980 double snorm(void)
1982 **********************************************************************
1985 (STANDARD-) N O R M A L DISTRIBUTION
1988 **********************************************************************
1989 **********************************************************************
1991 FOR DETAILS SEE:
1993 AHRENS, J.H. AND DIETER, U.
1994 EXTENSIONS OF FORSYTHE'S METHOD FOR RANDOM
1995 SAMPLING FROM THE NORMAL DISTRIBUTION.
1996 MATH. COMPUT., 27,124 (OCT. 1973), 927 - 937.
1998 ALL STATEMENT NUMBERS CORRESPOND TO THE STEPS OF ALGORITHM 'FL'
1999 (M=5) IN THE ABOVE PAPER (SLIGHTLY MODIFIED IMPLEMENTATION)
2001 Modified by Barry W. Brown, Feb 3, 1988 to use RANF instead of
2002 SUNIF. The argument IR thus goes away.
2004 **********************************************************************
2005 THE DEFINITIONS OF THE CONSTANTS A(K), D(K), T(K) AND
2006 H(K) ARE ACCORDING TO THE ABOVEMENTIONED ARTICLE
2009 static double a[32] = {
2010 0.0, 0.03917608550309, 0.07841241273311, 0.11776987457909,
2011 0.15731068461017, 0.19709908429430, 0.23720210932878, 0.27769043982157,
2012 0.31863936396437, 0.36012989178957, 0.40225006532172, 0.44509652498551,
2013 0.48877641111466, 0.53340970624127, 0.57913216225555, 0.62609901234641,
2014 0.67448975019607, 0.72451438349236, 0.77642176114792, 0.83051087820539,
2015 0.88714655901887, 0.94678175630104, 1.00999016924958, 1.07751556704027,
2016 1.15034938037600, 1.22985875921658, 1.31801089730353, 1.41779713799625,
2017 1.53412054435253, 1.67593972277344, 1.86273186742164, 2.15387469406144
2019 static double d[31] = {
2020 0.0, 0.0, 0.0, 0.0,
2021 0.0, 0.26368432217502, 0.24250845238097, 0.22556744380930,
2022 0.21163416577204, 0.19992426749317, 0.18991075842246, 0.18122518100691,
2023 0.17360140038056, 0.16684190866667, 0.16079672918053, 0.15534971747692,
2024 0.15040938382813, 0.14590257684509, 0.14177003276856, 0.13796317369537,
2025 0.13444176150074, 0.13117215026483, 0.12812596512583, 0.12527909006226,
2026 0.12261088288608, 0.12010355965651, 0.11774170701949, 0.11551189226063,
2027 0.11340234879117, 0.11140272044119, 0.10950385201710
2029 static double t[31] = {
2030 7.6738283767E-4, 2.30687039764E-3, 3.86061844387E-3, 5.43845406707E-3,
2031 7.05069876857E-3, 8.70839582019E-3, 1.042356984914E-2, 1.220953194966E-2,
2032 1.408124734637E-2, 1.605578804548E-2, 1.815290075142E-2, 2.039573175398E-2,
2033 2.281176732513E-2, 2.543407332319E-2, 2.830295595118E-2, 3.146822492920E-2,
2034 3.499233438388E-2, 3.895482964836E-2, 4.345878381672E-2, 4.864034918076E-2,
2035 5.468333844273E-2, 6.184222395816E-2, 7.047982761667E-2, 8.113194985866E-2,
2036 9.462443534514E-2, 0.11230007889456, 0.13649799954975, 0.17168856004707,
2037 0.22762405488269, 0.33049802776911, 0.58470309390507
2039 static double h[31] = {
2040 3.920617164634E-2, 3.932704963665E-2, 3.950999486086E-2, 3.975702679515E-2,
2041 4.007092772490E-2, 4.045532602655E-2, 4.091480886081E-2, 4.145507115859E-2,
2042 4.208311051344E-2, 4.280748137995E-2, 4.363862733472E-2, 4.458931789605E-2,
2043 4.567522779560E-2, 4.691571371696E-2, 4.833486978119E-2, 4.996298427702E-2,
2044 5.183858644724E-2, 5.401138183398E-2, 5.654656186515E-2, 5.953130423884E-2,
2045 6.308488965373E-2, 6.737503494905E-2, 7.264543556657E-2, 7.926471414968E-2,
2046 8.781922325338E-2, 9.930398323927E-2, 0.11555994154118, 0.14043438342816,
2047 0.18361418337460, 0.27900163464163, 0.70104742502766
2049 static long i;
2050 static double snorm,u,s,ustar,aa,w,y,tt;
2051 u = ranf();
2052 s = 0.0;
2053 if(u > 0.5) s = 1.0;
2054 u += (u-s);
2055 u = 32.0*u;
2056 i = (long) (u);
2057 if(i == 32) i = 31;
2058 if(i == 0) goto S100;
2060 START CENTER
2062 ustar = u-(double)i;
2063 aa = *(a+i-1);
2064 S40:
2065 if(ustar <= *(t+i-1)) goto S60;
2066 w = (ustar-*(t+i-1))**(h+i-1);
2067 S50:
2069 EXIT (BOTH CASES)
2071 y = aa+w;
2072 snorm = y;
2073 if(s == 1.0) snorm = -y;
2074 return snorm;
2075 S60:
2077 CENTER CONTINUED
2079 u = ranf();
2080 w = u*(*(a+i)-aa);
2081 tt = (0.5*w+aa)*w;
2082 goto S80;
2083 S70:
2084 tt = u;
2085 ustar = ranf();
2086 S80:
2087 if(ustar > tt) goto S50;
2088 u = ranf();
2089 if(ustar >= u) goto S70;
2090 ustar = ranf();
2091 goto S40;
2092 S100:
2094 START TAIL
2096 i = 6;
2097 aa = *(a+31);
2098 goto S120;
2099 S110:
2100 aa += *(d+i-1);
2101 i += 1;
2102 S120:
2103 u += u;
2104 if(u < 1.0) goto S110;
2105 u -= 1.0;
2106 S140:
2107 w = u**(d+i-1);
2108 tt = (0.5*w+aa)*w;
2109 goto S160;
2110 S150:
2111 tt = u;
2112 S160:
2113 ustar = ranf();
2114 if(ustar > tt) goto S50;
2115 u = ranf();
2116 if(ustar >= u) goto S150;
2117 u = ranf();
2118 goto S140;
2121 double fsign( double num, double sign )
2122 /* Transfers sign of argument sign to argument num */
2124 if ( ( sign>0.0f && num<0.0f ) || ( sign<0.0f && num>0.0f ) )
2125 return -num;
2126 else return num;
2129 /************************************************************************
2130 FTNSTOP:
2131 Prints msg to standard error and then exits
2132 ************************************************************************/
2133 void ftnstop(char* msg)
2134 /* msg - error message */
2136 if (msg != NULL) fprintf(stderr,"%s\n",msg);
2137 exit(0);