source: Sophya/trunk/SophyaLib/BaseTools/randinterf.cc@ 3735

Last change on this file since 3735 was 3615, checked in by cmv, 16 years ago

Modifs relatives a l'introduction de RandomGeneratorInterface + delete de srandgen.c, cmv 01/05/2009

File size: 30.0 KB
Line 
1#include "sopnamsp.h"
2#include "machdefs.h"
3#include <math.h>
4#include <stdlib.h>
5#include <sys/time.h>
6#include <time.h>
7#include <iostream>
8#include "pexceptions.h"
9
10#include "randinterf.h"
11
12namespace SOPHYA {
13
14//-------------------------------------------------------------------------------
15// ------ Definition d'interface des classes de generateurs de nombres aleatoires
16/*!
17 \class RandomGeneratorInterface
18 \ingroup BaseTools
19 \brief Base class for random number generators
20
21 This class defines the interface for random number generator classes and
22 implements the generation of some specific distributions (Gaussian, Poisson ...)
23 through generation of random number with a flat distribution in the range [0,1[.
24
25 The sub classes inheriting from this class should implement the Next() method.
26
27 This base class manages also a global instance of a default generator.
28
29 \sa frand01 drand01 frandpm1 drandpm1
30 \sa Gaussian Poisson
31
32*/
33
34
35RandomGeneratorInterface* RandomGeneratorInterface::gl_rndgen_p = NULL;
36
37/*!
38 \brief: static method to set or change the intance of the global Random Generator object
39
40 This method should be called during initialization, before any call to global
41 functions for random number generation. The rgp object should be created using new.
42*/
43void RandomGeneratorInterface::SetGlobalRandGenP(RandomGeneratorInterface* rgp)
44{
45 if (rgp == NULL) return;
46 if (gl_rndgen_p) delete gl_rndgen_p;
47 gl_rndgen_p = rgp;
48 return;
49}
50
51RandomGeneratorInterface::RandomGeneratorInterface()
52{
53 SelectGaussianAlgo();
54 SelectPoissonAlgo();
55 SelectExponentialAlgo();
56}
57
58
59RandomGeneratorInterface::~RandomGeneratorInterface(void)
60{
61 // rien a faire
62}
63
64void RandomGeneratorInterface::ShowRandom()
65{
66 cout<<"RandomGenerator is RandomGeneratorInterface i.e. UNDEFINED"<<endl;
67}
68
69/////////////////////////////////////////////////////////////////////////
70/////////////////////////////////////////////////////////////////////////
71/////////////////////////////////////////////////////////////////////////
72
73r_8 RandomGeneratorInterface::Next()
74{
75 printf("RandomGeneratorInterface::Next(): undefined code !!!\n");
76 throw MathExc("RandomGeneratorInterface::Next(): undefined code !!!");
77}
78
79/////////////////////////////////////////////////////////////////////////
80/////////////////////////////////////////////////////////////////////////
81/////////////////////////////////////////////////////////////////////////
82void RandomGeneratorInterface::GenerateSeedVector(int nseed,vector<uint_2>& seed,int lp)
83// renvoie un vecteur de nseed+2 entiers 32 bits
84// [0 - 2] = codage sur 48 bits du nombre (melange) de microsec depuis l'origine
85// [3 -> 3+ngene-1] = entiers aleatoires (poor man generator)
86//
87// L'initialiseur est donne par un codage du nombre de millisecondes
88// ecoulees depuis le 0 heure le 1er Janvier 1970 UTC (cf gettimeofday).
89// Seuls les 48 bits de poids faible sont retenus.
90// Un melange des bits est ensuite effectue pour que les 3 nombres
91// (unsigned short) d'initialisation ne soient pas trop semblables.
92// Le nombre le plus grand que l'on peut mettre
93// dans un entier unsigned de N bits est: 2^N-1
94// 48 bits -> 2^48-1 = 281474976710655 musec = 3257.8j = 8.9y
95// -> meme initialisation tous les 8.9 ans a 1 microsec pres !
96{
97 if(lp>0) cout<<"RandomGeneratorInterface::GenerateSeedVector: nseed="<<nseed<<endl;
98
99 // ---
100 // --- les deux premiers mots remplis avec le temps
101 // ---
102 // On recupere le temps ecoule depuis l'origine code en sec+musec
103 struct timeval now;
104 gettimeofday(&now,0);
105 // Calcul du temps ecoule depuis l'origine en microsecondes
106 uint_8 tmicro70 = (uint_8)now.tv_sec*(uint_8)1000000 + (uint_8)now.tv_usec;
107 if(lp>1) cout<<"."<<now.tv_sec<<" sec + "<<now.tv_usec<<" musec = "<<tmicro70<<" musec"<<endl;
108 // Remplissage du tableau de 48 bits
109 uint_2 b[48]; uint_8 tdum = tmicro70;
110 for(int ip=0;ip<48;ip++) {b[ip] = tdum&1; tdum = (tdum>>1);}
111 if(lp>2) {
112 cout<<"..b= ";
113 for(int ip=47;ip>=0;ip--) {cout<<b[ip]; if(ip%32==0 || ip%16==0) cout<<" ";}
114 cout<<endl;
115 }
116 // Melange des bits qui varient vite (poids faible, microsec)
117 // avec ceux variant lentement (poids fort, sec)
118 for(int ip=0;ip<16;ip++) {
119 if(ip%3==1) swap(b[ip],b[32+ip]);
120 else if(ip%3==2) swap(b[ip],b[16-ip]);
121 }
122 if(lp>2) {
123 cout<<"..b= ";
124 for(int ip=47;ip>=0;ip--) {cout<<b[ip]; if(ip%32==0 || ip%16==0) cout<<" ";}
125 cout<<endl;
126 }
127 // Remplissage
128 seed.resize(0);
129 for(int i=0;i<3;i++) {
130 seed.push_back(0);
131 uint_2 w = 1;
132 for(int ip=0;ip<16;ip++) {seed[i] += w*b[i*16+ip]; w *= 2;}
133 }
134 if(lp>0) cout<<"seed(time): "<<seed[0]<<" "<<seed[1]<<" "<<seed[2]<<endl;
135
136 // ---
137 // --- generation des nombres aleatoires complementaires (poor man generator)
138 // ---
139 //----------------------------------------------------------------------------//
140 // Ran088: L'Ecuyer's 1996 three-component Tausworthe generator "taus88"
141 // Returns an integer random number uniformly distributed within [0,4294967295]
142 // The period length is approximately 2^88 (which is 3*10^26).
143 // This generator is very fast and passes all standard statistical tests.
144 // Reference:
145 // (1) P. L'Ecuyer, Maximally equidistributed combined Tausworthe generators,
146 // Mathematics of Computation, 65, 203-213 (1996), see Figure 4.
147 // (2) recommended in:
148 // P. L'Ecuyer, Random number generation, chapter 4 of the
149 // Handbook on Simulation, Ed. Jerry Banks, Wiley, 1997.
150 //----------------------------------------------------------------------------//
151 if(nseed<=0) return;
152 // initialize seeds using the given seed value taking care of
153 // the requirements. The constants below are arbitrary otherwise
154 uint_4 seed0 = uint_4(tmicro70&0xFFFFFFFFULL);
155 if(lp>2) cout<<"seed0(time): "<<seed0<<endl;
156 uint_4 state_s1, state_s2, state_s3;
157 state_s1 = 1243598713U ^ seed0; if (state_s1 < 2) state_s1 = 1243598713U;
158 state_s2 = 3093459404U ^ seed0; if (state_s2 < 8) state_s2 = 3093459404U;
159 state_s3 = 1821928721U ^ seed0; if (state_s3 < 16) state_s3 = 1821928721U;
160 int nfill = 0, ico=0;
161 while(nfill<nseed) {
162 uint_4 s1 = state_s1, s2 = state_s2, s3 = state_s3;
163 // generate a random 32 bit number
164 s1 = ((s1 & -2) << 12) ^ (((s1 << 13) ^ s1) >> 19);
165 s2 = ((s2 & -8) << 4) ^ (((s2 << 2) ^ s2) >> 25);
166 s3 = ((s3 & -16) << 17) ^ (((s3 << 3) ^ s3) >> 11);
167 state_s1 = s1; state_s2 = s2; state_s3 = s3;
168 // le nombre aleatoire sur 32 bits est: s1^s2^s3
169 if(ico<15) {ico++; continue;} // des tirages blancs
170 uint_2 s = uint_2( (s1^s2^s3)&0xFFFFU );
171 seed.push_back(s);
172 if(lp>0) cout<<"seed(t88): "<<seed[3+nfill]<<endl;
173 nfill++;
174 }
175
176}
177
178void RandomGeneratorInterface::AutoInit(int lp)
179{
180 printf("RandomGeneratorInterface::AutoInit(): undefined code !!!\n");
181 throw MathExc("RandomGeneratorInterface::AutoInit(): undefined code !!!");
182}
183
184/////////////////////////////////////////////////////////////////////////
185/////////////////////////////////////////////////////////////////////////
186/////////////////////////////////////////////////////////////////////////
187
188r_8 RandomGeneratorInterface::Gaussian()
189{
190 switch (usegaussian_) {
191 case C_Gaussian_BoxMuller :
192 return GaussianBoxMuller();
193 break;
194 case C_Gaussian_RandLibSNorm :
195 return GaussianSNorm();
196 break;
197 case C_Gaussian_PolarBoxMuller :
198 return GaussianPolarBoxMuller();
199 break;
200 case C_Gaussian_RatioUnif :
201 return GaussianRatioUnif();
202 break;
203 case C_Gaussian_LevaRatioUnif :
204 return GaussianLevaRatioUnif();
205 break;
206 default:
207 return GaussianBoxMuller();
208 break;
209 }
210}
211
212//--- Generation de nombre aleatoires suivant une distribution gaussienne
213r_8 RandomGeneratorInterface::GaussianBoxMuller()
214{
215 r_8 A=Next();
216 while (A==0.) A=Next();
217 return sqrt(-2.*log(A))*cos(2.*M_PI*Next());
218}
219
220//-------------------------------------------
221// Adapte de ranlib float snorm()
222// http://orion.math.iastate.edu/burkardt/c_src/ranlib/ranlib.c
223/*
224**********************************************************************
225 (STANDARD-) N O R M A L DISTRIBUTION
226**********************************************************************
227
228 FOR DETAILS SEE:
229
230 AHRENS, J.H. AND DIETER, U.
231 EXTENSIONS OF FORSYTHE'S METHOD FOR RANDOM
232 SAMPLING FROM THE NORMAL DISTRIBUTION.
233 MATH. COMPUT., 27,124 (OCT. 1973), 927 - 937.
234
235 ALL STATEMENT NUMBERS CORRESPOND TO THE STEPS OF ALGORITHM 'FL'
236 (M=5) IN THE ABOVE PAPER (SLIGHTLY MODIFIED IMPLEMENTATION)
237
238 Modified by Barry W. Brown, Feb 3, 1988 to use RANF instead of
239 SUNIF. The argument IR thus goes away.
240
241**********************************************************************
242 THE DEFINITIONS OF THE CONSTANTS A(K), D(K), T(K) AND
243 H(K) ARE ACCORDING TO THE ABOVEMENTIONED ARTICLE
244*/
245static double a_snorm[32] = {
246 0.0,3.917609E-2,7.841241E-2,0.11777,0.1573107,0.1970991,0.2372021,0.2776904,
247 0.3186394,0.36013,0.4022501,0.4450965,0.4887764,0.5334097,0.5791322,
248 0.626099,0.6744898,0.7245144,0.7764218,0.8305109,0.8871466,0.9467818,
249 1.00999,1.077516,1.150349,1.229859,1.318011,1.417797,1.534121,1.67594,
250 1.862732,2.153875
251};
252static double d_snorm[31] = {
253 0.0,0.0,0.0,0.0,0.0,0.2636843,0.2425085,0.2255674,0.2116342,0.1999243,
254 0.1899108,0.1812252,0.1736014,0.1668419,0.1607967,0.1553497,0.1504094,
255 0.1459026,0.14177,0.1379632,0.1344418,0.1311722,0.128126,0.1252791,
256 0.1226109,0.1201036,0.1177417,0.1155119,0.1134023,0.1114027,0.1095039
257};
258static float t_snorm[31] = {
259 7.673828E-4,2.30687E-3,3.860618E-3,5.438454E-3,7.0507E-3,8.708396E-3,
260 1.042357E-2,1.220953E-2,1.408125E-2,1.605579E-2,1.81529E-2,2.039573E-2,
261 2.281177E-2,2.543407E-2,2.830296E-2,3.146822E-2,3.499233E-2,3.895483E-2,
262 4.345878E-2,4.864035E-2,5.468334E-2,6.184222E-2,7.047983E-2,8.113195E-2,
263 9.462444E-2,0.1123001,0.136498,0.1716886,0.2276241,0.330498,0.5847031
264};
265static float h_snorm[31] = {
266 3.920617E-2,3.932705E-2,3.951E-2,3.975703E-2,4.007093E-2,4.045533E-2,
267 4.091481E-2,4.145507E-2,4.208311E-2,4.280748E-2,4.363863E-2,4.458932E-2,
268 4.567523E-2,4.691571E-2,4.833487E-2,4.996298E-2,5.183859E-2,5.401138E-2,
269 5.654656E-2,5.95313E-2,6.308489E-2,6.737503E-2,7.264544E-2,7.926471E-2,
270 8.781922E-2,9.930398E-2,0.11556,0.1404344,0.1836142,0.2790016,0.7010474
271};
272r_8 RandomGeneratorInterface::GaussianSNorm()
273{
274long i;
275double snorm,u,s,ustar,aa,w,y,tt;
276 u = Next();
277 s = 0.0;
278 if(u > 0.5) s = 1.0;
279 u += (u-s);
280 u = 32.0*u;
281 i = (long) (u);
282 if(i == 32) i = 31;
283 if(i == 0) goto S100;
284/*
285 START CENTER
286*/
287 ustar = u-(double)i;
288 aa = *(a_snorm+i-1);
289S40:
290 if(ustar <= *(t_snorm+i-1)) goto S60;
291 w = (ustar-*(t_snorm+i-1))**(h_snorm+i-1);
292S50:
293/*
294 EXIT (BOTH CASES)
295*/
296 y = aa+w;
297 snorm = y;
298 if(s == 1.0) snorm = -y;
299 return snorm;
300S60:
301/*
302 CENTER CONTINUED
303*/
304 u = Next();
305 w = u*(*(a_snorm+i)-aa);
306 tt = (0.5*w+aa)*w;
307 goto S80;
308S70:
309 tt = u;
310 ustar = Next();
311S80:
312 if(ustar > tt) goto S50;
313 u = Next();
314 if(ustar >= u) goto S70;
315 ustar = Next();
316 goto S40;
317S100:
318/*
319 START TAIL
320*/
321 i = 6;
322 aa = *(a_snorm+31);
323 goto S120;
324S110:
325 aa += *(d_snorm+i-1);
326 i += 1;
327S120:
328 u += u;
329 if(u < 1.0) goto S110;
330 u -= 1.0;
331S140:
332 w = u**(d_snorm+i-1);
333 tt = (0.5*w+aa)*w;
334 goto S160;
335S150:
336 tt = u;
337S160:
338 ustar = Next();
339 if(ustar > tt) goto S50;
340 u = Next();
341 if(ustar >= u) goto S150;
342 u = Next();
343 goto S140;
344}
345
346r_8 RandomGeneratorInterface::GaussianPolarBoxMuller()
347{
348double x1,x2,w;
349do {
350 x1 = 2.0 * Next() - 1.0;
351 x2 = 2.0 * Next() - 1.0;
352 w = x1 * x1 + x2 * x2;
353 } while ( w >= 1.0 || w==0. );
354return x1 * sqrt(-2.0*log(w)/w);
355}
356
357static double s2se_RatioUnif=sqrt(2./M_E) , epm135_RatioUnif=exp(-1.35) , ep1q_RatioUnif=exp(1./4.);
358r_8 RandomGeneratorInterface::GaussianRatioUnif()
359{
360double u,v,x;
361while(true) {
362 do {u = Next();} while ( u == 0. );
363 v = (2.0*Next()-1.0)*s2se_RatioUnif;
364 x = v/u;
365 if(x*x <= 5.0-4.0*ep1q_RatioUnif*u) break;
366 if(x*x<4.0*epm135_RatioUnif/u+1.4)
367 if(v*v<-4.0*u*u*log(u)) break;
368}
369return x;
370}
371
372r_8 RandomGeneratorInterface::GaussianLevaRatioUnif()
373{
374double u,v,x,y,q;
375 do {
376 u = 1.-Next(); // in ]0,1]
377 v = Next()-0.5; // in [-0.5, 0.5[
378 v *= 1.7156;
379 x = u - 0.449871;
380 y = ((v<0)?-v:v) + 0.386595;
381 q = x*x + y*(0.19600*y - 0.25472*x);
382 } while( q>=0.27597 && (q>0.27846 || v*v>-4.0*u*u*log(u)) );
383 return v/u;
384}
385
386r_8 RandomGeneratorInterface::GaussianTail(double s)
387{
388 /* Returns a gaussian random variable larger than a
389 * This implementation does one-sided upper-tailed deviates.
390 */
391
392 if (s < 1)
393 {
394 /* For small s, use a direct rejection method. The limit s < 1
395 can be adjusted to optimise the overall efficiency */
396 double x;
397 do
398 {
399 x = Gaussian();
400 }
401 while (x < s);
402 return x;
403 }
404 else
405 {
406 /* Use the "supertail" deviates from the last two steps
407 * of Marsaglia's rectangle-wedge-tail method, as described
408 * in Knuth, v2, 3rd ed, pp 123-128. (See also exercise 11, p139,
409 * and the solution, p586.)
410 */
411 double u, v, x;
412 do
413 {
414 u = Next();
415 do
416 {
417 v = Next();
418 }
419 while (v == 0.0);
420 x = sqrt (s * s - 2 * log (v));
421 }
422 while (x * u > s);
423 return x;
424 }
425}
426
427/////////////////////////////////////////////////////////////////////////
428/////////////////////////////////////////////////////////////////////////
429/////////////////////////////////////////////////////////////////////////
430
431uint_8 RandomGeneratorInterface::Poisson(double mu, double mumax)
432{
433 switch (usepoisson_) {
434 case C_Poisson_Simple :
435 return PoissonSimple(mu,mumax);
436 break;
437 case C_Poisson_Ahrens :
438 return PoissonAhrens(mu);
439 break;
440 default:
441 return PoissonSimple(mu,mumax);
442 break;
443 }
444}
445
446
447//--- Generation de nombre aleatoires suivant une distribution de Poisson
448uint_8 RandomGeneratorInterface::PoissonSimple(double mu,double mumax)
449{
450 double pp,ppi,x;
451
452 if((mumax>0.)&&(mu>=mumax)) {
453 pp = sqrt(mu);
454 while( (x=pp*Gaussian()) < -mu );
455 return (uint_8)(mu+x+0.5);
456 }
457 else {
458 uint_8 n;
459 ppi = pp = exp(-mu);
460 x = Next();
461 n = 0;
462 while (x > ppi) {
463 n++;
464 pp = mu*pp/(double)n;
465 ppi += pp;
466 }
467 return n;
468 }
469 return 0; // pas necessaire ?
470}
471
472
473static double a0_poiahr = -0.5;
474static double a1_poiahr = 0.3333333;
475static double a2_poiahr = -0.2500068;
476static double a3_poiahr = 0.2000118;
477static double a4_poiahr = -0.1661269;
478static double a5_poiahr = 0.1421878;
479static double a6_poiahr = -0.1384794;
480static double a7_poiahr = 0.125006;
481static double fact_poiahr[10] = {
482 1.0,1.0,2.0,6.0,24.0,120.0,720.0,5040.0,40320.0,362880.0};
483uint_8 RandomGeneratorInterface::PoissonAhrens(double mu)
484/*
485**********************************************************************
486 long ignpoi(float mu)
487 GENerate POIsson random deviate
488 Function
489 Generates a single random deviate from a Poisson
490 distribution with mean AV.
491 Arguments
492 av --> The mean of the Poisson distribution from which
493 a random deviate is to be generated.
494 genexp <-- The random deviate.
495 Method
496 Renames KPOIS from TOMS as slightly modified by BWB to use RANF
497 instead of SUNIF.
498 For details see:
499 Ahrens, J.H. and Dieter, U.
500 Computer Generation of Poisson Deviates
501 From Modified Normal Distributions.
502 ACM Trans. Math. Software, 8, 2
503 (June 1982),163-179
504**********************************************************************
505**********************************************************************
506
507
508 P O I S S O N DISTRIBUTION
509
510
511**********************************************************************
512**********************************************************************
513
514 FOR DETAILS SEE:
515
516 AHRENS, J.H. AND DIETER, U.
517 COMPUTER GENERATION OF POISSON DEVIATES
518 FROM MODIFIED NORMAL DISTRIBUTIONS.
519 ACM TRANS. MATH. SOFTWARE, 8,2 (JUNE 1982), 163 - 179.
520
521 (SLIGHTLY MODIFIED VERSION OF THE PROGRAM IN THE ABOVE ARTICLE)
522
523**********************************************************************
524 INTEGER FUNCTION IGNPOI(IR,MU)
525 INPUT: IR=CURRENT STATE OF BASIC RANDOM NUMBER GENERATOR
526 MU=MEAN MU OF THE POISSON DISTRIBUTION
527 OUTPUT: IGNPOI=SAMPLE FROM THE POISSON-(MU)-DISTRIBUTION
528 MUPREV=PREVIOUS MU, MUOLD=MU AT LAST EXECUTION OF STEP P OR B.
529 TABLES: COEFFICIENTS A0-A7 FOR STEP F. FACTORIALS FACT
530 COEFFICIENTS A(K) - FOR PX = FK*V*V*SUM(A(K)*V**K)-DEL
531 SEPARATION OF CASES A AND B
532*/
533{
534uint_8 ignpoi,j,k,kflag,l,m;
535double b1,b2,c,c0,c1,c2,c3,d,del,difmuk,e,fk,fx,fy,g,omega,p,p0,px,py,q,s,
536 t,u,v,x,xx,pp[35];
537
538 if(mu < 10.0) goto S120;
539/*
540 C A S E A. (RECALCULATION OF S,D,L IF MU HAS CHANGED)
541*/
542 s = sqrt(mu);
543 d = 6.0*mu*mu;
544/*
545 THE POISSON PROBABILITIES PK EXCEED THE DISCRETE NORMAL
546 PROBABILITIES FK WHENEVER K >= M(MU). L=IFIX(MU-1.1484)
547 IS AN UPPER BOUND TO M(MU) FOR ALL MU >= 10 .
548*/
549 l = (uint_8) (mu-1.1484);
550/*
551 STEP N. NORMAL SAMPLE - SNORM(IR) FOR STANDARD NORMAL DEVIATE
552*/
553 g = mu+s*Gaussian();
554 if(g < 0.0) goto S20;
555 ignpoi = (uint_8) (g);
556/*
557 STEP I. IMMEDIATE ACCEPTANCE IF IGNPOI IS LARGE ENOUGH
558*/
559 if(ignpoi >= l) return ignpoi;
560/*
561 STEP S. SQUEEZE ACCEPTANCE - SUNIF(IR) FOR (0,1)-SAMPLE U
562*/
563 fk = (double)ignpoi;
564 difmuk = mu-fk;
565 u = Next();
566 if(d*u >= difmuk*difmuk*difmuk) return ignpoi;
567S20:
568/*
569 STEP P. PREPARATIONS FOR STEPS Q AND H.
570 (RECALCULATIONS OF PARAMETERS IF NECESSARY)
571 .3989423=(2*PI)**(-.5) .416667E-1=1./24. .1428571=1./7.
572 THE QUANTITIES B1, B2, C3, C2, C1, C0 ARE FOR THE HERMITE
573 APPROXIMATIONS TO THE DISCRETE NORMAL PROBABILITIES FK.
574 C=.1069/MU GUARANTEES MAJORIZATION BY THE 'HAT'-FUNCTION.
575*/
576 omega = 0.3989423/s;
577 b1 = 4.166667E-2/mu;
578 b2 = 0.3*b1*b1;
579 c3 = 0.1428571*b1*b2;
580 c2 = b2-15.0*c3;
581 c1 = b1-6.0*b2+45.0*c3;
582 c0 = 1.0-b1+3.0*b2-15.0*c3;
583 c = 0.1069/mu;
584 if(g < 0.0) goto S50;
585/*
586 'SUBROUTINE' F IS CALLED (KFLAG=0 FOR CORRECT RETURN)
587*/
588 kflag = 0;
589 goto S70;
590S40:
591/*
592 STEP Q. QUOTIENT ACCEPTANCE (RARE CASE)
593*/
594 if(fy-u*fy <= py*exp(px-fx)) return ignpoi;
595S50:
596/*
597 STEP E. EXPONENTIAL SAMPLE - SEXPO(IR) FOR STANDARD EXPONENTIAL
598 DEVIATE E AND SAMPLE T FROM THE LAPLACE 'HAT'
599 (IF T <= -.6744 THEN PK < FK FOR ALL MU >= 10.)
600*/
601 e = Exponential();
602 u = Next();
603 u += (u-1.0);
604 //t = 1.8+fsign(e,u);
605 t = 1.8 + (((u>0. && e<0.) || (u<0. && e>0.))?-e:e);
606 if(t <= -0.6744) goto S50;
607 ignpoi = (uint_8) (mu+s*t);
608 fk = (double)ignpoi;
609 difmuk = mu-fk;
610/*
611 'SUBROUTINE' F IS CALLED (KFLAG=1 FOR CORRECT RETURN)
612*/
613 kflag = 1;
614 goto S70;
615S60:
616/*
617 STEP H. HAT ACCEPTANCE (E IS REPEATED ON REJECTION)
618*/
619 if(c*fabs(u) > py*exp(px+e)-fy*exp(fx+e)) goto S50;
620 return ignpoi;
621S70:
622/*
623 STEP F. 'SUBROUTINE' F. CALCULATION OF PX,PY,FX,FY.
624 CASE IGNPOI .LT. 10 USES FACTORIALS FROM TABLE FACT
625*/
626 if(ignpoi >= 10) goto S80;
627 px = -mu;
628 py = pow(mu,(double)ignpoi)/ *(fact_poiahr+ignpoi);
629 goto S110;
630S80:
631/*
632 CASE IGNPOI .GE. 10 USES POLYNOMIAL APPROXIMATION
633 A0-A7 FOR ACCURACY WHEN ADVISABLE
634 .8333333E-1=1./12. .3989423=(2*PI)**(-.5)
635*/
636 del = 8.333333E-2/fk;
637 del -= (4.8*del*del*del);
638 v = difmuk/fk;
639 if(fabs(v) <= 0.25) goto S90;
640 px = fk*log(1.0+v)-difmuk-del;
641 goto S100;
642S90:
643 px = fk*v*v*(((((((a7_poiahr*v+a6_poiahr)*v+a5_poiahr)*v+a4_poiahr)*v+a3_poiahr)*v+a2_poiahr)*v+a1_poiahr)*v+a0_poiahr)-del;
644S100:
645 py = 0.3989423/sqrt(fk);
646S110:
647 x = (0.5-difmuk)/s;
648 xx = x*x;
649 fx = -0.5*xx;
650 fy = omega*(((c3*xx+c2)*xx+c1)*xx+c0);
651 if(kflag <= 0) goto S40;
652 goto S60;
653S120:
654/*
655 C A S E B. (START NEW TABLE AND CALCULATE P0 IF NECESSARY)
656*/
657// m = max(1L,(long) (mu));
658 m = (1ULL >= (uint_8)mu) ? 1ULL: (uint_8)mu;
659
660 l = 0;
661 p = exp(-mu);
662 q = p0 = p;
663S130:
664/*
665 STEP U. UNIFORM SAMPLE FOR INVERSION METHOD
666*/
667 u = Next();
668 ignpoi = 0;
669 if(u <= p0) return ignpoi;
670/*
671 STEP T. TABLE COMPARISON UNTIL THE END PP(L) OF THE
672 PP-TABLE OF CUMULATIVE POISSON PROBABILITIES
673 (0.458=PP(9) FOR MU=10)
674*/
675 if(l == 0) goto S150;
676 j = 1;
677//if(u > 0.458) j = min(l,m);
678 if(u > 0.458) j = ((l<=m)? l: m);
679 for(k=j; k<=l; k++) {
680 if(u <= *(pp+k-1)) goto S180;
681 }
682 if(l == 35) goto S130;
683S150:
684/*
685 STEP C. CREATION OF NEW POISSON PROBABILITIES P
686 AND THEIR CUMULATIVES Q=PP(K)
687*/
688 l += 1;
689 for(k=l; k<=35; k++) {
690 p = p*mu/(double)k;
691 q += p;
692 *(pp+k-1) = q;
693 if(u <= q) goto S170;
694 }
695 l = 35;
696 goto S130;
697S170:
698 l = k;
699S180:
700 ignpoi = k;
701 return ignpoi;
702}
703
704/////////////////////////////////////////////////////////////////////////
705/////////////////////////////////////////////////////////////////////////
706/////////////////////////////////////////////////////////////////////////
707
708r_8 RandomGeneratorInterface::Exponential()
709{
710 switch (useexpo_) {
711 case C_Exponential_Simple :
712 return ExpoSimple();
713 break;
714 case C_Exponential_Ahrens :
715 return ExpoAhrens();
716 break;
717 default:
718 return ExpoSimple();
719 break;
720 }
721}
722
723r_8 RandomGeneratorInterface::ExpoSimple(void)
724{
725 return -log(1.-Next());
726}
727
728
729static double q_expo[8] = {
730 0.6931472,0.9333737,0.9888778,0.9984959,0.9998293,0.9999833,0.9999986,1.0};
731r_8 RandomGeneratorInterface::ExpoAhrens(void)
732/*
733**********************************************************************
734**********************************************************************
735 (STANDARD-) E X P O N E N T I A L DISTRIBUTION
736**********************************************************************
737**********************************************************************
738
739 FOR DETAILS SEE:
740
741 AHRENS, J.H. AND DIETER, U.
742 COMPUTER METHODS FOR SAMPLING FROM THE
743 EXPONENTIAL AND NORMAL DISTRIBUTIONS.
744 COMM. ACM, 15,10 (OCT. 1972), 873 - 882.
745
746 ALL STATEMENT NUMBERS CORRESPOND TO THE STEPS OF ALGORITHM
747 'SA' IN THE ABOVE PAPER (SLIGHTLY MODIFIED IMPLEMENTATION)
748
749 Modified by Barry W. Brown, Feb 3, 1988 to use RANF instead of
750 SUNIF. The argument IR thus goes away.
751
752**********************************************************************
753 Q(N) = SUM(ALOG(2.0)**K/K!) K=1,..,N , THE HIGHEST N
754 (HERE 8) IS DETERMINED BY Q(N)=1.0 WITHIN STANDARD PRECISION
755*/
756{
757long i;
758double sexpo,a,u,ustar,umin;
759double *q1 = q_expo;
760 a = 0.0;
761 while((u=Next())==0.);
762 goto S30;
763S20:
764 a += *q1;
765S30:
766 u += u;
767 if(u <= 1.0) goto S20;
768 u -= 1.0;
769 if(u > *q1) goto S60;
770 sexpo = a+u;
771 return sexpo;
772S60:
773 i = 1;
774 ustar = Next();
775 umin = ustar;
776S70:
777 ustar = Next();
778 if(ustar < umin) umin = ustar;
779 i += 1;
780 if(u > *(q_expo+i-1)) goto S70;
781 sexpo = a+umin**q1;
782 return sexpo;
783}
784
785/////////////////////////////////////////////////////////////////////////
786/////////////////////////////////////////////////////////////////////////
787/////////////////////////////////////////////////////////////////////////
788
789int RandomGeneratorInterface::Gaussian2DRho(double &x,double &y,double mx,double my,double sx,double sy,double ro)
790/*
791++
792| Tirage de 2 nombres aleatoires x et y distribues sur une gaussienne 2D
793| de centre (mx,my), de coefficient de correlation rho (ro) et telle que
794| les sigmas finals des variables x et y soient sx,sy (ce sont
795| les valeurs des distributions marginales des variables aleatoires x et y
796| c'est a dire les sigmas des projections x et y de l'histogramme 2D
797| de la gaussienne). Retourne 0 si ok.
798|
799| - La densite de probabilite (normalisee a 1) sur laquelle on tire est:
800| N*exp[-0.5*{[(dx/sx)^2-2*ro/(sx*sy)*dx*dy+(dy/sy)^2]/(1-ro^2)}]
801| avec dx = x-mx, dy = y-my et N = 1/[2Pi*sx*sy*sqrt(1-ro^2)]
802| - Dans ce cas la distribution marginale est (ex en X):
803| 1/(sqrt(2Pi)*sx) * exp[-0.5*{dx^2/sx^2}]
804| - La matrice des covariances C des variables x,y est:
805| | sx^2 ro*sx*sy |
806| | | et det(C) = (1-ro^2)*sx^2*sy^2
807| | ro*sx*sy sy^2 |
808| - La matrice inverse C^(-1) est:
809| | 1/sx^2 -ro/(sx*sy) |
810| | | * 1/(1-ro^2)
811| | -ro/(sx*sy) 1/sy^2 |
812|
813| - Remarque:
814| le sigma que l'on obtient quand on fait une coupe de la gaussienne 2D
815| en y=0 (ou x=0) est: SX0(y=0) = sx*sqrt(1-ro^2) different de sx
816| SY0(x=0) = sy*sqrt(1-ro^2) different de sy
817| La distribution qui correspond a des sigmas SX0,SY0
818| pour les coupes en y=0,x=0 de la gaussienne 2D serait:
819| N*exp[-0.5*{ (dx/SX0)^2-2*ro/(SX0*SY0)*dx*dy+(dy/SY0)^2 }]
820| avec N = sqrt(1-ro^2)/(2Pi*SX0*SY0) et les variances
821| des variables x,y sont toujours
822| sx=SX0/sqrt(1-ro^2), sy=SY0/sqrt(1-ro^2)
823--
824*/
825{
826double a,b,sa;
827
828if( ro <= -1. || ro >= 1. ) return 1;
829
830while( (b=Flat01()) == 0. );
831b = sqrt(-2.*log(b));
832a = 2.*M_PI * Flat01();
833sa = sin(a);
834
835x = mx + sx*b*(sqrt(1.-ro*ro)*cos(a)+ro*sa);
836y = my + sy*b*sa;
837
838return 0;
839}
840
841void RandomGeneratorInterface::Gaussian2DAng(double &x,double &y,double mx,double my,double sa,double sb,double teta)
842/*
843++
844| Tirage de 2 nombres aleatoires x et y distribues sur une gaussienne 2D
845| de centre (x=mx,y=my), de sigmas grand axe et petit axe (sa,sb)
846| et dont le grand axe fait un angle teta (radian) avec l'axe des x.
847|
848| - La densite de probabilite (normalisee a 1) sur laquelle on tire est:
849| N*exp[-0.5*{ (A/sa)**2+(C/sc)**2 }], N=1/(2Pi*sa*sc)
850| ou A et B sont les coordonnees selon le grand axe et le petit axe
851| et teta = angle(x,A), le resultat subit ensuite une rotation d'angle teta.
852| - La matrice des covariances C des variables A,B est:
853| | sa^2 0 |
854| | | et det(C) = (1-ro^2)*sa^2*sb^2
855| | 0 sb^2 |
856| - La distribution x,y resultante est:
857| N*exp[-0.5*{[(dx/sx)^2-2*ro/(sx*sy)*dx*dy+(dy/sy)^2]/(1-ro^2)}]
858| ou N est donne dans NormCo et sx,sy,ro sont calcules a partir
859| de sa,sc,teta (voir fonctions paramga ou gaparam). La matrice des
860| covariances des variables x,y est donnee dans la fonction NormCo.
861--
862*/
863{
864double c,s,X,Y;
865
866while( (s = Flat01()) == 0. );
867s = sqrt(-2.*log(s));
868c = 2.*M_PI * Flat01();
869
870X = sa*s*cos(c);
871Y = sb*s*sin(c);
872
873c = cos(teta); s = sin(teta);
874x = mx + c*X - s*Y;
875y = my + s*X + c*Y;
876}
877
878} /* namespace SOPHYA */
879
880
881
882/////////////////////////////////////////////////////////////////
883/*
884**** Remarques sur complex< r_8 > ComplexGaussian(double sig) ****
885
886--- variables gaussiennes x,y independantes
887x gaussien: pdf f(x) = 1/(sqrt(2Pi) Sx) exp(-(x-Mx)^2/(2 Sx^2))
888y gaussien: pdf f(y) = 1/(sqrt(2Pi) Sy) exp(-(y-My)^2/(2 Sy^2))
889x,y independants --> pdf f(x,y) = f(x) f(y)
890On a:
891 <x> = Integrate[x*f(x)] = Mx
892 <x^2> = Integrate[x^2*f(x)] = Mx^2 + Sx^2
893
894--- On cherche la pdf g(r,t) du module et de la phase
895 x = r cos(t) , y = r sin(t)
896 r=sqrt(x^2+y^2 , t=atan2(y,x)
897 (r,t) --> (x,y): le Jacobien = r
898
899 g(r,t) = r f(x,y) = r f(x) f(y)
900 = r/(2Pi Sx Sy) exp(-(x-Mx)^2/(2 Sx^2)) exp(-(y-My)^2/(2 Sy^2))
901
902- Le cas general est complique
903 (cf D.Pelat cours DEA "bruits et signaux" section 4.5)
904
905- Cas ou "Mx = My = 0" et "Sx = Sy = S"
906 c'est la pdf du module et de la phase d'un nombre complexe
907 dont les parties reelles et imaginaires sont independantes
908 et sont distribuees selon des gaussiennes de variance S^2
909 g(r,t) = r/(2Pi S^2) exp(-r^2/(2 S^2))
910 La distribution de "r" est donc:
911 g(r) = Integrate[g(r,t),{t,0,2Pi}]
912 = r/S^2 exp(-r^2/(2 S^2))
913 La distribution de "t" est donc:
914 g(t) = Integrate[g(r,t),{r,0,Infinity}]
915 = 1 / 2Pi (distribution uniforme sur [0,2Pi[)
916 Les variables aleatoires r,t sont independantes:
917 g(r,t) = g(r) g(t)
918On a:
919 <r> = Integrate[r*g(r)] = sqrt(PI/2)*S
920 <r^2> = Integrate[r^2*g(r)] = 2*S^2
921 <r^3> = Integrate[r^3*g(r)] = 3*sqrt(PI/2)*S^3
922 <r^4> = Integrate[r^4*g(r)] = 8*S^4
923
924- Attention:
925La variable complexe "c = x+iy = r*exp(i*t)" ainsi definie verifie:
926 <|c|^2> = <c c*> = <x^2+y^2> = <r^2> = 2 S^2
927Si on veut generer une variable complexe gaussienne telle que
928 <c c*> = s^2 alors il faut prendre S = s/sqrt(2) comme argument
929
930*/
Note: See TracBrowser for help on using the repository browser.