Main Page   Namespace List   Class Hierarchy   Alphabetical List   Compound List   File List   Namespace Members   Compound Members  

trnglib.cc

00001 // Time-stamp: <Donnerstag, 18.04.2002, 11:48:43; edited by heiko>
00002 // 
00003 // Tina's random number generators TRNG
00004 //
00005 // useful functions for Tina's random number generators 
00006 //
00007 // Copyright (C) 2001, 2002 Heiko Bauke
00008 //
00009 // heiko.bauke@student.uni-magdeburg.de
00010 //
00011 // TRNG is free software; you can redistribute it and/or
00012 // modify it under the terms of the GNU General Public License
00013 // as published by the Free Software Foundation. This program
00014 // is distributed WITHOUT ANY WARRANTY; without even the implied
00015 // warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
00016 // See the GNU General Public License for more details.
00017 //
00018 // ---------------------------------------------------------------------
00019 
00020 #include <trng.h>
00021 #include <trnglib.h>
00022 #include <cmath>
00023 
00025 
00030 const char * TRNG::version(void) {
00031   return "TRNG V"TRNG_VERSION"."TRNG_REVISION;
00032 }
00033 
00034 // ---------------------------------------------------------------------
00035 
00037 
00046 long TRNG::modulo_invers(long a, long m) {
00047   if (a<=0l || m<=1)
00048     throw error("invalid argument in TNRG::modulo_invers(long, long)");
00049   long long al=static_cast<long long>(a);
00050   long e=m-2l;
00051   long long p=1ll;
00052   while (e>0l) {
00053     if ((e&1l)==1l) {
00054       p*=al;
00055       p%=m;
00056     }
00057     al*=al;
00058     al%=m;
00059     e>>=1;
00060   }
00061   return static_cast<long>(p);
00062 }
00063 
00064 // ---------------------------------------------------------------------
00065   
00067 
00096 void TRNG::gauss(std::vector<long> &a, std::vector<long> &b, long m) {
00097   if (a.size()!=b.size()*b.size() || a.size()==0 || b.size()==0) 
00098     throw TRNG::error("wrong matrix size in TRNG::gauss");
00099   // initialize indices
00100   int n=b.size();
00101   std::vector<long> p(n);
00102   for (int i=0; i<n; ++i)
00103     p[i]=i;
00104   // make matrix diagonal
00105   for (int i=0; i<n-1; ++i) {
00106     if (a[n*p[i]+i]==0l) {
00107       // swap
00108       int j=i+1;
00109       while (j<n && a[n*p[j]+i]==0l) 
00110         j++;
00111       if (j==n)
00112         throw TRNG::error("singular matrix in TRNG::gauss");
00113       long t=p[i];  p[i]=p[j];  p[j]=t;
00114     }
00115     long t=TRNG::modulo_invers(a[n*p[i]+i], m);
00116     for (int j=i; j<n; ++j) 
00117       a[n*p[i]+j]=static_cast<long>
00118         ((static_cast<long long>(a[n*p[i]+j])*
00119           static_cast<long long>(t))%m);
00120     b[p[i]]=static_cast<long>
00121       ((static_cast<long long>(b[p[i]])*
00122         static_cast<long long>(t))%m);
00123     for (int j=i+1; j<n; ++j) {
00124       if (a[n*p[j]+i]!=0l) {
00125         t=TRNG::modulo_invers(a[n*p[j]+i], m);
00126         for (int k=i; k<n; ++k) {
00127           a[n*p[j]+k]=
00128             static_cast<long>
00129             ((static_cast<long long>(a[n*p[j]+k])*
00130               static_cast<long long>(t))%m);
00131           a[n*p[j]+k]-=a[n*p[i]+k];
00132           if (a[n*p[j]+k]<0l)
00133             a[n*p[j]+k]+=m;
00134         }
00135         b[p[j]]=static_cast<long>
00136           ((static_cast<long long>(b[p[j]])*
00137             static_cast<long long>(t))%m);
00138         b[p[j]]-=b[p[i]];
00139         if (b[p[j]]<0l)
00140           b[p[j]]+=m;
00141       }
00142     }
00143   }
00144   // solve diagonal system
00145   long t=TRNG::modulo_invers(a[n*p[n-1]+n-1], m);
00146   a[n*p[n-1]+n-1]=1l;
00147   b[p[n-1]]=static_cast<long>
00148     ((static_cast<long long>(b[p[n-1]])*
00149       static_cast<long long>(t))%m);
00150   for (int i=n-2; i>=0; --i)
00151     for (int j=i+1; j<n; ++j) {
00152       b[p[i]]-=static_cast<long>
00153         ((static_cast<long long>(a[n*p[i]+j])*
00154           static_cast<long long>(b[p[j]]))%m);
00155       if (b[p[i]]<0l)
00156         b[p[i]]+=m;
00157     }
00158   // sort
00159   for (int i=0; i<n; ++i) 
00160     p[i]=b[p[i]];
00161   for (int i=0; i<n; ++i) 
00162     b[i]=p[i];
00163 }
00164 
00165 // ---------------------------------------------------------------------
00166 
00168 
00179 void TRNG::matrix_mult(const std::vector<long> &a, const std::vector<long> &b,
00180                        std::vector<long> &c, long m) {
00181   if (a.size()!=b.size())
00182     throw TRNG::error("different sized matrices in TRNG::matrix_mult");
00183   int n=static_cast<int>(std::sqrt(static_cast<double>(b.size())));
00184   for (int i=0; i<n; ++i)
00185     for (int j=0; j<n; ++j) {
00186       long long t=0ll;
00187       for (int k=0; k<n; ++k) {
00188         t+=(static_cast<long long>(a[j*n+k])*
00189             static_cast<long long>(b[k*n+i]))%m;
00190         if (t>=m)
00191           t-=m;
00192       }
00193       c[j*n+i]=static_cast<long>(t);
00194     }
00195 }
00196 
00197 // ---------------------------------------------------------------------
00198 
00200 
00212 void TRNG::matrix_vec_mult(const std::vector<long> &a, 
00213                            const std::vector<long> &b, 
00214                            std::vector<long> &c, long m) {
00215   if (a.size()!=b.size()*b.size())
00216     throw TRNG::error("different sized vectors in TRNG::matrix_vec_mult");
00217   int n=b.size();
00218   for (int j=0; j<n; ++j) {
00219     long long t=0ll;
00220     for (int k=0; k<n; ++k) {
00221       t+=(static_cast<long long>(a[j*n+k])*
00222           static_cast<long long>(b[k]))%m;
00223       if (t>=m)
00224         t-=m;
00225     }
00226     c[j]=static_cast<long>(t);
00227   }
00228 }
00229 
00230 // ---------------------------------------------------------------------
00231 
00232 // for special functions calculation see
00233 // Samuel S. M. Wong Computational Methods in Physics and Engineering
00234 
00235 // ---------------------------------------------------------------------
00236 
00238 
00248 double TRNG::ln_Gamma(double x) {
00249   if (x<=0)
00250     throw TRNG::error("invalid (nonpositive) argument for TRNG::lnGamma");
00251   const double gamma=7.0;
00252   const int N=9;
00253   const double ln_sqrt_2_pi=.9189385332046725; // ln(sqrt(2 Pi))
00254   double t1=x-0.5;
00255   double t2=t1+gamma;
00256   /* coefficients for gamma=7, kmax=8  Lanczos method */
00257   const double lanczos_7_c[N]={0.99999999999980993227684700473478,
00258                                676.520368121885098567009190444019,
00259                                -1259.13921672240287047156078755283,
00260                                771.3234287776530788486528258894,
00261                                -176.61502916214059906584551354,
00262                                12.507343278686904814458936853,
00263                                -0.13857109526572011689554707,
00264                                9.984369578019570859563e-6,
00265                                1.50563273514931155834e-7};
00266   t1=ln_sqrt_2_pi+t1*std::log(t2)-t2;
00267   t2=lanczos_7_c[0];
00268   for (int i=1; i<N; ++i)
00269     t2+=lanczos_7_c[i]/x++;
00270   return t1+std::log(t2);
00271 }
00272 
00273 // ---------------------------------------------------------------------
00274 
00276 
00285 double TRNG::Gamma(double x) {
00286   if (x<=0.0)
00287     throw TRNG::error("invalid (nonpositive) argument for TRNG::Gamma");
00288   // use a interpolating polynom for 1<=x<2
00289   if (x<15.0) {
00290     double f, z;
00291     // Gamma(x)=Gamma(x+1)/x
00292     if (x<1.0) {
00293       f=1.0/x;
00294       z=x;
00295     }
00296     if (1.0<=x && x<2.0) {
00297       z=x-1.0;
00298       f=1.0;
00299     }
00300     // Gamma(x+1)=x*Gamma(x)
00301     if (x>=2.0) {
00302       f=1.0;
00303       z=x-1.0;
00304       do {
00305         f=f*z;
00306         z--;
00307       } while (z>=1.0);
00308     }
00309     return f*(1.0+
00310               (-0.5772152952998961+
00311                (0.9890418621994901+
00312                 (-0.9072503205615771+
00313                  (0.9796139097861202+
00314                   (-0.9693649942410396+
00315                    (0.9409502806615237+
00316                     (-0.8401327394958836+
00317                      (0.6506941079199484+
00318                       (-0.4028899587819087+
00319                        (0.1812400864733248+
00320                         (-0.5156307405048575E-1
00321                          +0.6876141229E-2
00322                          *z)*z)*z)*z)*z)*z)*z)*z)*z)*z)*z)*z);
00323   }
00324   // Stirling's approximation
00325   double t=x*x;
00326   return std::exp(-x)*std::pow(x, x-0.5)*2.5066282746310*
00327     (1.0
00328      +1.0/12.0/x
00329      +1/288.0/t
00330      -139.0/51840.0/(t*x)
00331      -571.0/2488320.0/(t*t));
00332 }
00333 
00334 // ---------------------------------------------------------------------
00335 
00337 
00348 double TRNG::Gamma_P(double a, double x) {
00349   if (x<0.0 || a<=0.0)
00350     throw TRNG::error("invalid (nonpositive) argument for TRNG::GammaP");
00351   double t;
00352   t=TRNG::Gamma(a);
00353   if (x<a+1.0)
00354     return TRNG::Gamma_ser(a, x)/t;
00355   else
00356     return 1.0-TRNG::Gamma_cf(a, x)/t;
00357 }
00358 
00359 // ---------------------------------------------------------------------
00360 
00362 
00373 double TRNG::Gamma_Q(double a, double x) {
00374   if (x<0.0 || a<=0.0)
00375     throw TRNG::error("invalid (nonpositive) argument for TRNG::GammaQ");
00376   double t;
00377   t=TRNG::Gamma(a);
00378   if (x<a+1.0)
00379     return 1.0-TRNG::Gamma_ser(a, x)/t;
00380   else
00381     return TRNG::Gamma_cf(a, x)/t;
00382 }
00383 
00384 // ---------------------------------------------------------------------
00385 
00387 
00398 double TRNG::incomp_Gamma(double a, double x) {
00399   if (x<0.0 || a<=0.0)
00400     throw TRNG::error("invalid (nonpositive) argument for TRNG::incomp_Gamma");
00401   if (x<a+1.0)
00402     return TRNG::Gamma_ser(a, x);
00403   else
00404     return TRNG::Gamma(a)-TRNG::Gamma_cf(a, x);
00405 }
00406 
00407 // ---------------------------------------------------------------------
00408 
00410 
00421 double TRNG::comp_incomp_Gamma(double a, double x) {
00422   if (x<0.0 || a<=0.0)
00423     throw TRNG::error("invalid (nonpositive) argument for TRNG::comp_incomp_Gamma");
00424   if (x<a+1.0)
00425     return TRNG::Gamma(a)-TRNG::Gamma_ser(a, x);
00426   else
00427     return TRNG::Gamma_cf(a, x);
00428 }
00429 
00430 // ---------------------------------------------------------------------
00431 
00433 
00442 double TRNG::Gamma_ser(double a, double x) {
00443   const int itmax=150;
00444   const double eps=1e-12;
00445   int i;
00446   double sum, xx, n; 
00447   if (x<0.0)
00448     throw TRNG::error("x less than 0 in routine TRNG::Gamma_ser");
00449   if (a<=0.0)
00450     throw TRNG::error("a less or equal than 0 in routine TRNG::Gamma_ser");
00451   if (x==0.0)
00452     return 0.0;
00453   xx=1.0/a;
00454   n=a;
00455   sum=xx;
00456   i=0;
00457   do {
00458     ++n;
00459     ++i;
00460     xx*=x/n;
00461     sum+=xx;
00462   } while (std::fabs(xx)>eps*std::fabs(sum) && i<itmax);
00463   if (i==itmax)
00464     throw TRNG::error("convergence problem in TRNG::Gamma_ser");
00465   return std::exp(-x+a*std::log(x))*sum;
00466 }
00467 
00468 // ---------------------------------------------------------------------
00469 
00471 
00480 double TRNG::Gamma_cf(double a, double x) {
00481   // complementary incomplete Gamma function's continued fraction  
00482   // representation for x>a+1
00483   const double itmax=75.0;
00484   const double eps=1.0e-12;
00485   const double min=1.0e-45;
00486   double ai, bi, ci, di, del, h, i; 
00487   if (x<0.0)
00488     throw TRNG::error("x less than 0 in routine TRNG::Gamma_cf");
00489   if (a<=0.0)
00490     throw TRNG::error("a less or equal than 0 in routine TRNG::Gamma_cf");
00491   // Set up for evaluating continued fraction by modied Lentz's method
00492   bi=x+1.0-a; 
00493   ci=1.0/min; 
00494   di=1.0/bi; 
00495   h=di;
00496   i=0.0;
00497   do { 
00498     // Iterate to convergence. 
00499     ++i;
00500     ai=-i*(i-a); 
00501     bi+=2.0;
00502     di=ai*di+bi;
00503     if (std::fabs(di)<min) 
00504       di=min; 
00505     ci=bi+ai/ci; 
00506     if (std::fabs(ci)<min)
00507       ci=min;
00508     di=1.0/di;
00509     del=di*ci;
00510     h*=del;
00511   } while ((std::fabs(del-1.0)>eps) && i<itmax);
00512   if (i==itmax)
00513     throw TRNG::error("a too large or convergence problem in TRNG::Gamma_cf"); 
00514   return std::exp(-x+a*std::log(x))*h;
00515 }
00516 
00517 // ---------------------------------------------------------------------
00518 
00520 
00527 double TRNG::ln_factorial(long n) {
00528   if (n<0l)
00529     throw TRNG::error("parameter less than zero in TRNG::ln_factorial");
00530   const long n_ln_fac_tab=32l;
00531   static double ln_fac_tab[n_ln_fac_tab]={0.0, 0.0};
00532   static long in_tab_count=1l;
00533   if (n>=n_ln_fac_tab)
00534     return TRNG::ln_Gamma(static_cast<double>(n)+1.0);
00535   if (n>in_tab_count) {
00536     for (long i=in_tab_count+1; i<=n; ++i)
00537       ln_fac_tab[i]=ln_fac_tab[i-1l]+std::log(static_cast<double>(i));
00538     in_tab_count=n;
00539   }
00540   return ln_fac_tab[n];
00541 }
00542 
00543 // ---------------------------------------------------------------------
00544 
00546 
00553 long TRNG::binomial_coeff(long n, long k) {
00554   if (n<0 || k<0 || k>n) 
00555     return 0l;
00556   return static_cast<long>
00557     (std::exp(TRNG::ln_factorial(n)-TRNG::ln_factorial(k)-
00558               TRNG::ln_factorial(n-k))+0.5);
00559 }
00560 
00561 // ---------------------------------------------------------------------
00562 
00564 
00573 double TRNG::errf(double x) {
00574   return (x>0.0) ? 
00575     0.5*TRNG::Gamma_P(0.5, 0.5*x*x) : 
00576     -0.5*TRNG::Gamma_P(0.5, 0.5*x*x);
00577 }
00578 
00579 // ---------------------------------------------------------------------
00580 
00582 
00592 double TRNG::chi_square_test(const std::vector<double> &prob, 
00593                              const std::vector<double> &observ) {
00594   double chi2, n, t1, t2;
00595   if (prob.size()!=observ.size())
00596     throw TRNG::error("different sized vectors in TRNG::chi_square_test");
00597   n=0.0;
00598   chi2=0.0;
00599   for (unsigned int i=0; i<prob.size(); ++i) {
00600     n+=observ[i];
00601     if (observ[i]<5.0)
00602       throw TRNG::error("not enough observations in TRNG::chi_square_test");
00603   }
00604   for (unsigned int i=0; i<prob.size(); ++i) {
00605     t1=prob[i]*n;
00606     t2=observ[i]-t1;
00607     chi2+=t2*t2/t1;
00608   }
00609   return chi2;
00610 }
00611 
00612 // ---------------------------------------------------------------------
00613 
00615 
00624 double TRNG::chi_square_prob(double chi2, long df) {
00625   if (df<1l)
00626     throw TRNG::error("degrees of freedom less than one in TRNG::chi_square_prob");
00627   double ddf;
00628   ddf=static_cast<double>(df);
00629   if (chi2==0.0)
00630     return 1.0;
00631   if (df<256l)
00632     return TRNG::Gamma_Q(0.5*ddf, 0.5*chi2);
00633   // see Knuth, The Art of Computer Programming II, page 41
00634   return 0.5-TRNG::errf(-3.0/4.0*std::sqrt(2.0)*std::sqrt(ddf)+
00635                         1.0/4.0*std::sqrt(-6.0*ddf+24.0*chi2+16.0));
00636 }
00637 
00638 // ---------------------------------------------------------------------
00639 
00641 
00651 double TRNG::Stirling_num2(long n, long m) {
00652   // see Knuth, The Art of Computer Programming I
00653   if (n<0l || m<0l || m>n)
00654     return 0.0;
00655   if (m==n)
00656     return 1.0;
00657   if (m==0l)
00658     return 0.0;
00659   if (m==1l)
00660     return 1.0;
00661   if (m==2l)
00662     return static_cast<double>((1ll<<(n-1l))-1ll);
00663   return m*TRNG::Stirling_num2(n-1l, m)+TRNG::Stirling_num2(n-1l, m-1l);
00664 }
00665 
00666 // ---------------------------------------------------------------------
00667 
00669 
00691 double TRNG::Student_t(double p, long nu, bool symmetric) {
00692   const double pi=3.14159265358979324;
00693   const double sqrt2=1.41421356237309505;
00694   const double two_div_sqrt_Pi=1.12837916709551257;
00695   double t, dt, t2, g1, g2, g3, g4, m;
00696   int i;
00697   bool less_0_5;
00698   if (nu<1l)
00699     throw TRNG::error("less than one degree of freedom in TRNG::Student_t");
00700   if (p<=0.0 || p>=1.0)
00701     throw TRNG::error("probability out of range in TRNG::Student_t");
00702   m=static_cast<double>(nu);
00703   if (symmetric)
00704     p+=0.5*(1.0-p);
00705   // exact formulas from Statistik. Lehr- und Handbuch der angewandten 
00706   // Statistik., Joachim Hartung, Bärbel Elpelt, Karl-Heinz Klösener 
00707   // R. Oldenbourg Verlag, 1998, page 892
00708   if (nu==1l) {
00709     t=std::tan(pi*(p-0.5));
00710   }
00711   if (nu==2l) {
00712     t=(2.0*p-1.0);
00713     t*=sqrt2/std::sqrt(1.0-t*t); 
00714   }
00715   if (nu>2l) {
00716     // use symmetry t(1-p, m)=-t(p, m)
00717     less_0_5=(p<0.5);
00718     if (less_0_5)
00719       p=1.0-p;
00720     // claculate t(p, m) for m=oo with Newton's method 
00721     t=0.0;
00722     i=0;
00723     do {
00724       dt=-(2.0*TRNG::errf(sqrt2*t)+1.0-2.0*p)/(two_div_sqrt_Pi*std::exp(-t*t));
00725       t+=dt;
00726     } while (std::fabs(dt/t)>1e-14 && ++i<12);
00727     // use approximation form Handbook of Mathematical Functions, 
00728     // With Formulas, Graphs, and Mathematical Tables, Milton Abramowitz 
00729     // Dover Publications, Inc., page 949
00730     t*=sqrt2;
00731     t2=t*t;
00732     g1=1.0/4.0*(1.0+t2)*t;
00733     g2=1.0/96.0*(3.0+(16.0+5.0*t2)*t2)*t;
00734     g3=1.0/384.0*(-15.0+(17.0+.0*(19.0+3.0*t2)*t2)*t2)*t;
00735     g4=1.0/92160.0*(-945.0+(-1920.0+(1482.0+(776.0+79*t2)*t2)*t2)*t2)*t;
00736     t+=g1/m+g2/m/m+g3/m/m/m+g4/m/m/m/m;
00737     if (less_0_5)
00738       t*=-1.0;
00739   }
00740   return t;
00741 }
00742 
00743 // ---------------------------------------------------------------------
00744 
00746 
00762 long TRNG::find_interval(const std::vector<double> &borders, 
00763                          const double x) {
00764   long num_classes=borders.size()+1;
00765   if (num_classes==1 || x<=borders[0])
00766     return 0;
00767   if (borders[num_classes-2]<x)
00768     return num_classes-1;
00769   long i1=0;
00770   long i2=num_classes-2;
00771   while (i2-i1>1) {
00772     long i3=(i2+i1)/2;
00773     if (x<=borders[i3])
00774       i2=i3;
00775     else
00776       i1=i3;
00777   }
00778   return i2;
00779 }
00780 
00781 // ---------------------------------------------------------------------
00782 
00784 
00790 double TRNG::uniform_pdf(double x) {
00791   return uniformco_pdf(x);
00792 }
00793 
00795 
00801 double TRNG::uniform_pdf(double x, double a, double b) {
00802   return uniformco_pdf(x, a, b);
00803 }
00804 
00805 // ---------------------------------------------------------------------
00806 
00808 
00814 double TRNG::uniformco_pdf(double x) {
00815   if (0.0<=x && x<1.0)
00816     return 1.0;
00817   else
00818     return 0.0;
00819 }
00820 
00822 
00828 double TRNG::uniformco_pdf(double x, double a, double b) {
00829   if (a<=x && x<b)
00830     return 1.0/(b-a);
00831   else
00832     return 0.0;
00833 }
00834 
00835 // ---------------------------------------------------------------------
00836 
00838 
00844 double TRNG::uniformcc_pdf(double x) {
00845   if (0.0<=x && x<=1.0)
00846     return 1.0;
00847   else
00848     return 0.0;
00849 }
00850 
00852 
00859 double TRNG::uniformcc_pdf(double x, double a, double b) {
00860   if (a<=x && x<=b)
00861     return 1.0/(b-a);
00862   else
00863     return 0.0;
00864 }
00865 
00866 // ---------------------------------------------------------------------
00867 
00869 
00875 double TRNG::uniformoc_pdf(double x) {
00876   if (0.0<x && x<=1.0)
00877     return 1.0;
00878   else
00879     return 0.0;
00880 }
00881 
00883 
00889 double TRNG::uniformoc_pdf(double x, double a, double b) {
00890   if (a<x && x<=b)
00891     return 1.0/(b-a);
00892   else
00893     return 0.0;
00894 }
00895 
00896 // ---------------------------------------------------------------------
00897 
00899 
00905 double TRNG::uniformoo_pdf(double x) {
00906   if (0.0<x && x<1.0)
00907     return 1.0;
00908   else
00909     return 0.0;
00910 }
00911 
00913 
00919 double TRNG::uniformoo_pdf(double x, double a, double b) {
00920   if (a<x && x<b)
00921     return 1.0/(b-a);
00922   else
00923     return 0.0;
00924 }
00925 
00926 // ---------------------------------------------------------------------
00927 
00929 
00944 double TRNG::normal_dist_pdf(double x, double sigma=1.0, double mu=0.0) {
00945   const double one_sqrt_2_pi=.3989422804014327; // 1/sqrt(2*Pi)
00946   if (sigma<=0.0)
00947     throw TRNG::error("negative or zero standard deviation in TRNG::normal_dist_pdf");
00948   return one_sqrt_2_pi/sigma*std::exp(-(x-mu)*(x-mu)/(2.0*sigma*sigma));
00949 }
00950 
00951 // ---------------------------------------------------------------------
00952 
00954 
00971 double TRNG::exp_dist_pdf(double x, double mu=1.0) {
00972   if (mu<=0.0)
00973     throw TRNG::error("negative or zero parameter in TRNG::exp_dist_pdf");
00974   if (x>=0.0)
00975     return std::exp(-x/mu)/mu;
00976   else
00977     return 0.0;
00978 }
00979 
00980 // ---------------------------------------------------------------------
00981 
00983 
00996 double TRNG::laplace_dist_pdf(double x, double a=1.0) {
00997   if (a<=0.0)
00998     throw TRNG::error("negative or zero parameter in TRNG::exp_dist_pdf");
00999   return 0.5/a*std::exp(-std::fabs(x)/a);
01000 }
01001 
01002 // ---------------------------------------------------------------------
01003 
01005 
01022 double TRNG::tent_dist_pdf(double x, double a=1.0) {
01023   if (a<=0.0)
01024     throw TRNG::error("negative or zero parameter in TRNG::exp_dist_pdf");
01025   if (-a<x && x<=0.0)
01026     return (x+a)/(a*a);
01027   if (0.0<x && x<=a)
01028     return (a-x)/(a*a);
01029   return 0.0;
01030 }
01031 
01032 // ---------------------------------------------------------------------
01033 
01035 
01053 double TRNG::Gamma_dist_pdf(double x, double a, double b) {
01054   if (a<=0.0 || b<=0.0)
01055     throw TRNG::error("parameter less than one in TRNG::Gamma_dist_pdf");
01056   if (x<0.0) 
01057     return 0.0;
01058   if (x==0.0)
01059     if (a==1.0)
01060       return 1.0/b ;
01061     else
01062       return 0.0;
01063   if (a==1.0) 
01064     return std::exp(-x/b)/b ;
01065   return std::pow(x, a-1.0)*std::exp(-x/b)/TRNG::Gamma(a)/std::pow(b, a);
01066 }
01067 
01068 // ---------------------------------------------------------------------
01069 
01071 
01090 double TRNG::Beta_dist_pdf(double x, double a, double b) {
01091   if (a<=0.0 || b<=0.0)
01092     throw TRNG::error("parameter less than one in TRNG::Gamma_dist");
01093   if (x<0.0 || x>1.0)
01094     return 0.0;
01095   return std::exp(TRNG::ln_Gamma(a+b)-TRNG::ln_Gamma(a)-TRNG::ln_Gamma(b))*
01096     std::pow(x, a-1.0)*std::pow(1.0-x, b-1.0);
01097 }
01098 
01099 // ---------------------------------------------------------------------
01100 
01102 
01122 double TRNG::chi_square_dist_pdf(double x, double nu) {
01123   if (nu<1.0)
01124     throw TRNG::error("parameter less than one in TRNG::chi_square_dist_pdf");
01125   if (x<=0.0)
01126     return 0.0;
01127   return std::exp((0.5*nu-1.0)*std::log(0.5*x)-
01128                   0.5*x-TRNG::ln_Gamma(0.5*nu))*0.5;
01129 }
01130 
01131 // ---------------------------------------------------------------------
01132 
01134 
01149 double TRNG::Student_t_dist_pdf(double x, double nu) {
01150   if (nu<=0.0)
01151     throw TRNG::error("parameter less than or equal zero in TRNG::Student_t_dist_pdf");
01152   const double pi=3.14159265358979324;
01153   double t=(nu+1.0)*0.5;
01154   return std::exp(TRNG::ln_Gamma(t)-TRNG::ln_Gamma(nu/2.0))/
01155     std::sqrt(pi*nu)*std::pow((1.0+x*x/nu), -t);
01156 }
01157 
01158 // ---------------------------------------------------------------------
01159 
01161 
01174 double TRNG::binomial_dist_pdf(long k, long n, double p=0.5) {
01175   if (p<=0.0 || p>1.0)
01176     throw TRNG::error("probability <=0.0 or >1.0 in TRNG::binomial_dist_pdf");
01177   if (n<1l)
01178     throw TRNG::error("less than one trail TRNG::binomial_dist_pdf");
01179   if (k<0 || k>n)
01180     return 0.0;
01181   return TRNG::binomial_coeff(n, k)*
01182     std::pow(p, static_cast<double>(k))*
01183     std::pow(1.0-p, static_cast<double>(n-k));
01184 }
01185 
01186 // ---------------------------------------------------------------------
01187 
01189 
01207 double TRNG::poisson_dist_pdf(long k, double mu) {
01208   if (k<0l)
01209     return 0.0;
01210   return std::pow(mu, static_cast<double>(k))*
01211     std::exp(-mu-TRNG::ln_factorial(k));
01212 }
01213 
01214 // ---------------------------------------------------------------------
01215 
01217 
01229 double TRNG::geometric_dist_pdf(long k, double q) {
01230   if (q<=0.0 ||q>1.0)
01231     throw TRNG::error("parameter out of range in TRNG::RNG::geometric_dist_pdf");
01232   if (k<=0l)
01233     return 0.0;
01234   if (k==1l)
01235     return q;
01236   return q*std::pow(1.0-q, static_cast<double>(k)-1.0);
01237 }

Generated at Tue Apr 30 12:03:09 2002 for Tina's Random Number Generators by doxygen1.2.3 written by Dimitri van Heesch, © 1997-2000