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

trnglib.cc

00001 // Time-stamp: <Samstag, 25.01.2003, 14:36:02; 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@physik.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");
00049   long temp, q, flast=0, f=1, m1=m;
00050   while (a>1) {
00051     temp=m1%a; 
00052     q=m1/a;
00053     m1=a;  a=temp;  temp=f; 
00054     f=flast-q*f; 
00055     flast=temp;
00056   }
00057   if (a==0) 
00058     throw error("no inversive in TNRG::modulo_invers");
00059   return f>=0 ? f : f+m;
00060 }
00061 
00062 // ---------------------------------------------------------------------
00063   
00065 
00094 void TRNG::gauss(std::vector<long> &a, std::vector<long> &b, long m) {
00095   if (a.size()!=b.size()*b.size() || a.size()==0 || b.size()==0) 
00096     throw TRNG::error("wrong matrix size in TRNG::gauss");
00097   // initialize indices
00098   int n=b.size();
00099   int rank=0;
00100   std::vector<long> p(n);
00101   for (int i=0; i<n; ++i)
00102     p[i]=i;
00103   // make matrix triangular
00104   for (int i=0; i<n; ++i) {
00105     // search for a pivot element
00106     if (a[n*p[i]+i]==0l) {
00107       // swap rows
00108       int j=i+1;
00109       while (j<n && a[n*p[j]+i]==0l) 
00110         j++;
00111       if (j<n) {
00112         long t=p[i];  p[i]=p[j];  p[j]=t;
00113       }
00114     }
00115     // is rank small?
00116     if (a[n*p[i]+i]==0l)
00117       break;
00118     ++rank;
00119     long t=TRNG::modulo_invers(a[n*p[i]+i], m);
00120     for (int j=i; j<n; ++j) 
00121       a[n*p[i]+j]=static_cast<long>
00122         ((static_cast<long long>(a[n*p[i]+j])*
00123           static_cast<long long>(t))%m);
00124     b[p[i]]=static_cast<long>
00125       ((static_cast<long long>(b[p[i]])*
00126         static_cast<long long>(t))%m);
00127     for (int j=i+1; j<n; ++j) {
00128       if (a[n*p[j]+i]!=0l) {
00129         t=TRNG::modulo_invers(a[n*p[j]+i], m);
00130         for (int k=i; k<n; ++k) {
00131           a[n*p[j]+k]=
00132             static_cast<long>
00133             ((static_cast<long long>(a[n*p[j]+k])*
00134               static_cast<long long>(t))%m);
00135           a[n*p[j]+k]-=a[n*p[i]+k];
00136           if (a[n*p[j]+k]<0l)
00137             a[n*p[j]+k]+=m;
00138         }
00139         b[p[j]]=static_cast<long>
00140           ((static_cast<long long>(b[p[j]])*
00141             static_cast<long long>(t))%m);
00142         b[p[j]]-=b[p[i]];
00143         if (b[p[j]]<0l)
00144           b[p[j]]+=m;
00145       }
00146     }
00147   }
00148   // test if a solution exists
00149   for (int i=rank; i<n; ++i) 
00150     if (b[p[i]]!=0l)
00151       throw TRNG::error("equations system has no solution TRNG::gauss");
00152   // solve triangular system
00153   for (int i=n-2; i>=0; --i)
00154     for (int j=i+1; j<n; ++j) {
00155       b[p[i]]-=static_cast<long>
00156         ((static_cast<long long>(a[n*p[i]+j])*
00157           static_cast<long long>(b[p[j]]))%m);
00158       if (b[p[i]]<0l)
00159         b[p[i]]+=m;
00160     }
00161   // sort
00162   for (int i=0; i<n; ++i) 
00163     p[i]=b[p[i]];
00164   for (int i=0; i<n; ++i) 
00165     b[i]=p[i];
00166 }
00167 
00168 
00169 // ---------------------------------------------------------------------
00170 
00172 
00183 void TRNG::matrix_mult(const std::vector<long> &a, const std::vector<long> &b,
00184                        std::vector<long> &c, long m) {
00185   if (a.size()!=b.size())
00186     throw TRNG::error("different sized matrices in TRNG::matrix_mult");
00187   int n=static_cast<int>(std::sqrt(static_cast<double>(b.size())));
00188   c.resize(n*n);
00189   for (int i=0; i<n; ++i)
00190     for (int j=0; j<n; ++j) {
00191       long long t=0ll;
00192       for (int k=0; k<n; ++k) {
00193         t+=(static_cast<long long>(a[j*n+k])*
00194             static_cast<long long>(b[k*n+i]))%m;
00195         if (t>=m)
00196           t-=m;
00197       }
00198       c[j*n+i]=static_cast<long>(t);
00199     }
00200 }
00201 
00202 // ---------------------------------------------------------------------
00203 
00205 
00217 void TRNG::matrix_vec_mult(const std::vector<long> &a, 
00218                            const std::vector<long> &b, 
00219                            std::vector<long> &c, long m) {
00220   if (a.size()!=b.size()*b.size())
00221     throw TRNG::error("different sized vectors in TRNG::matrix_vec_mult");
00222   int n=b.size();
00223   c.resize(n);
00224   for (int j=0; j<n; ++j) {
00225     long long t=0ll;
00226     for (int k=0; k<n; ++k) {
00227       t+=(static_cast<long long>(a[j*n+k])*
00228           static_cast<long long>(b[k]))%m;
00229       if (t>=m)
00230         t-=m;
00231     }
00232     c[j]=static_cast<long>(t);
00233   }
00234 }
00235 
00236 // ---------------------------------------------------------------------
00237 
00238 // for special functions calculation see
00239 // Samuel S. M. Wong Computational Methods in Physics and Engineering
00240 
00241 // ---------------------------------------------------------------------
00242 
00244 
00254 double TRNG::ln_Gamma(double x) {
00255   if (x<=0)
00256     throw TRNG::error("invalid (nonpositive) argument for TRNG::lnGamma");
00257   const double gamma=7.0;
00258   const int N=9;
00259   const double ln_sqrt_2_pi=.9189385332046725; // ln(sqrt(2 Pi))
00260   double t1=x-0.5;
00261   double t2=t1+gamma;
00262   /* coefficients for gamma=7, kmax=8  Lanczos method */
00263   const double lanczos_7_c[N]={0.99999999999980993227684700473478,
00264                                676.520368121885098567009190444019,
00265                                -1259.13921672240287047156078755283,
00266                                771.3234287776530788486528258894,
00267                                -176.61502916214059906584551354,
00268                                12.507343278686904814458936853,
00269                                -0.13857109526572011689554707,
00270                                9.984369578019570859563e-6,
00271                                1.50563273514931155834e-7};
00272   t1=ln_sqrt_2_pi+t1*std::log(t2)-t2;
00273   t2=lanczos_7_c[0];
00274   for (int i=1; i<N; ++i)
00275     t2+=lanczos_7_c[i]/x++;
00276   return t1+std::log(t2);
00277 }
00278 
00279 // ---------------------------------------------------------------------
00280 
00282 
00291 double TRNG::Gamma(double x) {
00292   if (x<=0.0)
00293     throw TRNG::error("invalid (nonpositive) argument for TRNG::Gamma");
00294   // use a interpolating polynom for 1<=x<2
00295   if (x<15.0) {
00296     double f, z;
00297     // Gamma(x)=Gamma(x+1)/x
00298     if (x<1.0) {
00299       f=1.0/x;
00300       z=x;
00301     }
00302     if (1.0<=x && x<2.0) {
00303       z=x-1.0;
00304       f=1.0;
00305     }
00306     // Gamma(x+1)=x*Gamma(x)
00307     if (x>=2.0) {
00308       f=1.0;
00309       z=x-1.0;
00310       do {
00311         f=f*z;
00312         z--;
00313       } while (z>=1.0);
00314     }
00315     return f*(1.0+
00316               (-0.5772152952998961+
00317                (0.9890418621994901+
00318                 (-0.9072503205615771+
00319                  (0.9796139097861202+
00320                   (-0.9693649942410396+
00321                    (0.9409502806615237+
00322                     (-0.8401327394958836+
00323                      (0.6506941079199484+
00324                       (-0.4028899587819087+
00325                        (0.1812400864733248+
00326                         (-0.5156307405048575E-1
00327                          +0.6876141229E-2
00328                          *z)*z)*z)*z)*z)*z)*z)*z)*z)*z)*z)*z);
00329   }
00330   // Stirling's approximation
00331   double t=x*x;
00332   return std::exp(-x)*std::pow(x, x-0.5)*2.5066282746310*
00333     (1.0
00334      +1.0/12.0/x
00335      +1/288.0/t
00336      -139.0/51840.0/(t*x)
00337      -571.0/2488320.0/(t*t));
00338 }
00339 
00340 // ---------------------------------------------------------------------
00341 
00343 
00354 double TRNG::Gamma_P(double a, double x) {
00355   if (x<0.0 || a<=0.0)
00356     throw TRNG::error("invalid (nonpositive) argument for TRNG::GammaP");
00357   double t;
00358   t=TRNG::Gamma(a);
00359   if (x<a+1.0)
00360     return TRNG::Gamma_ser(a, x)/t;
00361   else
00362     return 1.0-TRNG::Gamma_cf(a, x)/t;
00363 }
00364 
00365 // ---------------------------------------------------------------------
00366 
00368 
00379 double TRNG::Gamma_Q(double a, double x) {
00380   if (x<0.0 || a<=0.0)
00381     throw TRNG::error("invalid (nonpositive) argument for TRNG::GammaQ");
00382   double t;
00383   t=TRNG::Gamma(a);
00384   if (x<a+1.0)
00385     return 1.0-TRNG::Gamma_ser(a, x)/t;
00386   else
00387     return TRNG::Gamma_cf(a, x)/t;
00388 }
00389 
00390 // ---------------------------------------------------------------------
00391 
00393 
00404 double TRNG::incomp_Gamma(double a, double x) {
00405   if (x<0.0 || a<=0.0)
00406     throw TRNG::error("invalid (nonpositive) argument for TRNG::incomp_Gamma");
00407   if (x<a+1.0)
00408     return TRNG::Gamma_ser(a, x);
00409   else
00410     return TRNG::Gamma(a)-TRNG::Gamma_cf(a, x);
00411 }
00412 
00413 // ---------------------------------------------------------------------
00414 
00416 
00427 double TRNG::comp_incomp_Gamma(double a, double x) {
00428   if (x<0.0 || a<=0.0)
00429     throw TRNG::error("invalid (nonpositive) argument for TRNG::comp_incomp_Gamma");
00430   if (x<a+1.0)
00431     return TRNG::Gamma(a)-TRNG::Gamma_ser(a, x);
00432   else
00433     return TRNG::Gamma_cf(a, x);
00434 }
00435 
00436 // ---------------------------------------------------------------------
00437 
00439 
00448 double TRNG::Gamma_ser(double a, double x) {
00449   const int itmax=150;
00450   const double eps=1e-12;
00451   int i;
00452   double sum, xx, n; 
00453   if (x<0.0)
00454     throw TRNG::error("x less than 0 in routine TRNG::Gamma_ser");
00455   if (a<=0.0)
00456     throw TRNG::error("a less or equal than 0 in routine TRNG::Gamma_ser");
00457   if (x==0.0)
00458     return 0.0;
00459   xx=1.0/a;
00460   n=a;
00461   sum=xx;
00462   i=0;
00463   do {
00464     ++n;
00465     ++i;
00466     xx*=x/n;
00467     sum+=xx;
00468   } while (std::fabs(xx)>eps*std::fabs(sum) && i<itmax);
00469   if (i==itmax)
00470     throw TRNG::error("convergence problem in TRNG::Gamma_ser");
00471   return std::exp(-x+a*std::log(x))*sum;
00472 }
00473 
00474 // ---------------------------------------------------------------------
00475 
00477 
00486 double TRNG::Gamma_cf(double a, double x) {
00487   // complementary incomplete Gamma function's continued fraction  
00488   // representation for x>a+1
00489   const double itmax=75.0;
00490   const double eps=1.0e-12;
00491   const double min=1.0e-45;
00492   double ai, bi, ci, di, del, h, i; 
00493   if (x<0.0)
00494     throw TRNG::error("x less than 0 in routine TRNG::Gamma_cf");
00495   if (a<=0.0)
00496     throw TRNG::error("a less or equal than 0 in routine TRNG::Gamma_cf");
00497   // Set up for evaluating continued fraction by modied Lentz's method
00498   bi=x+1.0-a; 
00499   ci=1.0/min; 
00500   di=1.0/bi; 
00501   h=di;
00502   i=0.0;
00503   do { 
00504     // Iterate to convergence. 
00505     ++i;
00506     ai=-i*(i-a); 
00507     bi+=2.0;
00508     di=ai*di+bi;
00509     if (std::fabs(di)<min) 
00510       di=min; 
00511     ci=bi+ai/ci; 
00512     if (std::fabs(ci)<min)
00513       ci=min;
00514     di=1.0/di;
00515     del=di*ci;
00516     h*=del;
00517   } while ((std::fabs(del-1.0)>eps) && i<itmax);
00518   if (i==itmax)
00519     throw TRNG::error("a too large or convergence problem in TRNG::Gamma_cf"); 
00520   return std::exp(-x+a*std::log(x))*h;
00521 }
00522 
00523 // ---------------------------------------------------------------------
00524 
00526 
00533 double TRNG::ln_factorial(long n) {
00534   if (n<0l)
00535     throw TRNG::error("parameter less than zero in TRNG::ln_factorial");
00536   const long n_ln_fac_tab=32l;
00537   static double ln_fac_tab[n_ln_fac_tab]={0.0, 0.0};
00538   static long in_tab_count=1l;
00539   if (n>=n_ln_fac_tab)
00540     return TRNG::ln_Gamma(static_cast<double>(n)+1.0);
00541   if (n>in_tab_count) {
00542     for (long i=in_tab_count+1; i<=n; ++i)
00543       ln_fac_tab[i]=ln_fac_tab[i-1l]+std::log(static_cast<double>(i));
00544     in_tab_count=n;
00545   }
00546   return ln_fac_tab[n];
00547 }
00548 
00549 // ---------------------------------------------------------------------
00550 
00552 
00559 long TRNG::binomial_coeff(long n, long k) {
00560   if (n<0 || k<0 || k>n) 
00561     return 0l;
00562   return static_cast<long>
00563     (std::exp(TRNG::ln_factorial(n)-TRNG::ln_factorial(k)-
00564               TRNG::ln_factorial(n-k))+0.5);
00565 }
00566 
00567 // ---------------------------------------------------------------------
00568 
00570 
00579 double TRNG::errf(double x) {
00580   return (x>0.0) ? 
00581     0.5*TRNG::Gamma_P(0.5, 0.5*x*x) : 
00582     -0.5*TRNG::Gamma_P(0.5, 0.5*x*x);
00583 }
00584 
00585 // ---------------------------------------------------------------------
00586 
00588 
00598 double TRNG::chi_square_test(const std::vector<double> &prob, 
00599                              const std::vector<double> &observ) {
00600   double chi2, n, t1, t2;
00601   if (prob.size()!=observ.size())
00602     throw TRNG::error("different sized vectors in TRNG::chi_square_test");
00603   n=0.0;
00604   chi2=0.0;
00605   for (unsigned int i=0; i<prob.size(); ++i) {
00606     n+=observ[i];
00607     if (observ[i]<5.0)
00608       throw TRNG::error("not enough observations in TRNG::chi_square_test");
00609   }
00610   for (unsigned int i=0; i<prob.size(); ++i) {
00611     t1=prob[i]*n;
00612     t2=observ[i]-t1;
00613     chi2+=t2*t2/t1;
00614   }
00615   return chi2;
00616 }
00617 
00618 // ---------------------------------------------------------------------
00619 
00621 
00630 double TRNG::chi_square_prob(double chi2, long df) {
00631   if (df<1l)
00632     throw TRNG::error("degrees of freedom less than one in TRNG::chi_square_prob");
00633   double ddf;
00634   ddf=static_cast<double>(df);
00635   if (chi2==0.0)
00636     return 1.0;
00637   if (df<256l)
00638     return TRNG::Gamma_Q(0.5*ddf, 0.5*chi2);
00639   // see Knuth, The Art of Computer Programming II, page 41
00640   return 0.5-TRNG::errf(-3.0/4.0*std::sqrt(2.0)*std::sqrt(ddf)+
00641                         1.0/4.0*std::sqrt(-6.0*ddf+24.0*chi2+16.0));
00642 }
00643 
00644 // ---------------------------------------------------------------------
00645 
00647 
00657 double TRNG::Stirling_num2(long n, long m) {
00658   // see Knuth, The Art of Computer Programming I
00659   if (n<0l || m<0l || m>n)
00660     return 0.0;
00661   if (m==n)
00662     return 1.0;
00663   if (m==0l)
00664     return 0.0;
00665   if (m==1l)
00666     return 1.0;
00667   if (m==2l)
00668     return static_cast<double>((1ll<<(n-1l))-1ll);
00669   return m*TRNG::Stirling_num2(n-1l, m)+TRNG::Stirling_num2(n-1l, m-1l);
00670 }
00671 
00672 // ---------------------------------------------------------------------
00673 
00675 
00697 double TRNG::Student_t(double p, long nu, bool symmetric) {
00698   const double pi=3.14159265358979324;
00699   const double sqrt2=1.41421356237309505;
00700   const double two_div_sqrt_Pi=1.12837916709551257;
00701   double t, dt, t2, g1, g2, g3, g4, m;
00702   int i;
00703   bool less_0_5;
00704   if (nu<1l)
00705     throw TRNG::error("less than one degree of freedom in TRNG::Student_t");
00706   if (p<=0.0 || p>=1.0)
00707     throw TRNG::error("probability out of range in TRNG::Student_t");
00708   m=static_cast<double>(nu);
00709   if (symmetric)
00710     p+=0.5*(1.0-p);
00711   // exact formulas from Statistik. Lehr- und Handbuch der angewandten 
00712   // Statistik., Joachim Hartung, Bärbel Elpelt, Karl-Heinz Klösener 
00713   // R. Oldenbourg Verlag, 1998, page 892
00714   if (nu==1l) {
00715     t=std::tan(pi*(p-0.5));
00716   }
00717   if (nu==2l) {
00718     t=(2.0*p-1.0);
00719     t*=sqrt2/std::sqrt(1.0-t*t); 
00720   }
00721   if (nu>2l) {
00722     // use symmetry t(1-p, m)=-t(p, m)
00723     less_0_5=(p<0.5);
00724     if (less_0_5)
00725       p=1.0-p;
00726     // claculate t(p, m) for m=oo with Newton's method 
00727     t=0.0;
00728     i=0;
00729     do {
00730       dt=-(2.0*TRNG::errf(sqrt2*t)+1.0-2.0*p)/(two_div_sqrt_Pi*std::exp(-t*t));
00731       t+=dt;
00732     } while (std::fabs(dt/t)>1e-14 && ++i<12);
00733     // use approximation form Handbook of Mathematical Functions, 
00734     // With Formulas, Graphs, and Mathematical Tables, Milton Abramowitz 
00735     // Dover Publications, Inc., page 949
00736     t*=sqrt2;
00737     t2=t*t;
00738     g1=1.0/4.0*(1.0+t2)*t;
00739     g2=1.0/96.0*(3.0+(16.0+5.0*t2)*t2)*t;
00740     g3=1.0/384.0*(-15.0+(17.0+.0*(19.0+3.0*t2)*t2)*t2)*t;
00741     g4=1.0/92160.0*(-945.0+(-1920.0+(1482.0+(776.0+79*t2)*t2)*t2)*t2)*t;
00742     t+=g1/m+g2/m/m+g3/m/m/m+g4/m/m/m/m;
00743     if (less_0_5)
00744       t*=-1.0;
00745   }
00746   return t;
00747 }
00748 
00749 // ---------------------------------------------------------------------
00750 
00752 
00768 long TRNG::find_interval(const std::vector<double> &borders, const double x) {
00769   long num_classes=borders.size()+1;
00770   if (num_classes==1 || x<=borders[0])
00771     return 0;
00772   if (borders[num_classes-2]<x)
00773     return num_classes-1;
00774   long i1=0;
00775   long i2=num_classes-2;
00776   while (i2-i1>1) {
00777     long i3=(i2+i1)/2;
00778     if (x<=borders[i3])
00779       i2=i3;
00780     else
00781       i1=i3;
00782   }
00783   return i2;
00784 }
00785 
00786 // ---------------------------------------------------------------------
00787 
00789 
00795 double TRNG::uniform_pdf(double x) {
00796   return uniformco_pdf(x);
00797 }
00798 
00800 
00806 double TRNG::uniform_pdf(double x, double a, double b) {
00807   return uniformco_pdf(x, a, b);
00808 }
00809 
00810 // ---------------------------------------------------------------------
00811 
00813 
00819 double TRNG::uniformco_pdf(double x) {
00820   if (0.0<=x && x<1.0)
00821     return 1.0;
00822   else
00823     return 0.0;
00824 }
00825 
00827 
00833 double TRNG::uniformco_pdf(double x, double a, double b) {
00834   if (a<=x && x<b)
00835     return 1.0/(b-a);
00836   else
00837     return 0.0;
00838 }
00839 
00840 // ---------------------------------------------------------------------
00841 
00843 
00849 double TRNG::uniformcc_pdf(double x) {
00850   if (0.0<=x && x<=1.0)
00851     return 1.0;
00852   else
00853     return 0.0;
00854 }
00855 
00857 
00864 double TRNG::uniformcc_pdf(double x, double a, double b) {
00865   if (a<=x && x<=b)
00866     return 1.0/(b-a);
00867   else
00868     return 0.0;
00869 }
00870 
00871 // ---------------------------------------------------------------------
00872 
00874 
00880 double TRNG::uniformoc_pdf(double x) {
00881   if (0.0<x && x<=1.0)
00882     return 1.0;
00883   else
00884     return 0.0;
00885 }
00886 
00888 
00894 double TRNG::uniformoc_pdf(double x, double a, double b) {
00895   if (a<x && x<=b)
00896     return 1.0/(b-a);
00897   else
00898     return 0.0;
00899 }
00900 
00901 // ---------------------------------------------------------------------
00902 
00904 
00910 double TRNG::uniformoo_pdf(double x) {
00911   if (0.0<x && x<1.0)
00912     return 1.0;
00913   else
00914     return 0.0;
00915 }
00916 
00918 
00924 double TRNG::uniformoo_pdf(double x, double a, double b) {
00925   if (a<x && x<b)
00926     return 1.0/(b-a);
00927   else
00928     return 0.0;
00929 }
00930 
00931 // ---------------------------------------------------------------------
00932 
00934 
00949 double TRNG::normal_dist_pdf(double x, double sigma=1.0, double mu=0.0) {
00950   const double one_sqrt_2_pi=.3989422804014327; // 1/sqrt(2*Pi)
00951   if (sigma<=0.0)
00952     throw TRNG::error("negative or zero standard deviation in TRNG::normal_dist_pdf");
00953   return one_sqrt_2_pi/sigma*std::exp(-(x-mu)*(x-mu)/(2.0*sigma*sigma));
00954 }
00955 
00956 // ---------------------------------------------------------------------
00957 
00959 
00977 double TRNG::exp_dist_pdf(double x, double mu=1.0) {
00978   if (mu<=0.0)
00979     throw TRNG::error("negative or zero parameter in TRNG::exp_dist_pdf");
00980   if (x>=0.0)
00981     return std::exp(-x/mu)/mu;
00982   else
00983     return 0.0;
00984 }
00985 
00986 // ---------------------------------------------------------------------
00987 
00989 
01003 double TRNG::laplace_dist_pdf(double x, double a=1.0) {
01004   if (a<=0.0)
01005     throw TRNG::error("negative or zero parameter in TRNG::exp_dist_pdf");
01006   return 0.5/a*std::exp(-std::fabs(x)/a);
01007 }
01008 
01009 // ---------------------------------------------------------------------
01010 
01012 
01030 double TRNG::tent_dist_pdf(double x, double a=1.0) {
01031   if (a<=0.0)
01032     throw TRNG::error("negative or zero parameter in TRNG::exp_dist_pdf");
01033   if (-a<x && x<=0.0)
01034     return (x+a)/(a*a);
01035   if (0.0<x && x<=a)
01036     return (a-x)/(a*a);
01037   return 0.0;
01038 }
01039 
01040 // ---------------------------------------------------------------------
01041 
01043 
01061 double TRNG::Gamma_dist_pdf(double x, double a, double b) {
01062   if (a<=0.0 || b<=0.0)
01063     throw TRNG::error("parameter less than one in TRNG::Gamma_dist_pdf");
01064   if (x<0.0) 
01065     return 0.0;
01066   if (x==0.0)
01067     if (a==1.0)
01068       return 1.0/b ;
01069     else
01070       return 0.0;
01071   if (a==1.0) 
01072     return std::exp(-x/b)/b ;
01073   return std::pow(x, a-1.0)*std::exp(-x/b)/TRNG::Gamma(a)/std::pow(b, a);
01074 }
01075 
01076 // ---------------------------------------------------------------------
01077 
01079 
01098 double TRNG::Beta_dist_pdf(double x, double a, double b) {
01099   if (a<=0.0 || b<=0.0)
01100     throw TRNG::error("parameter less than one in TRNG::Gamma_dist");
01101   if (x<0.0 || x>1.0)
01102     return 0.0;
01103   return std::exp(TRNG::ln_Gamma(a+b)-TRNG::ln_Gamma(a)-TRNG::ln_Gamma(b))*
01104     std::pow(x, a-1.0)*std::pow(1.0-x, b-1.0);
01105 }
01106 
01107 // ---------------------------------------------------------------------
01108 
01110 
01130 double TRNG::chi_square_dist_pdf(double x, double nu) {
01131   if (nu<1.0)
01132     throw TRNG::error("parameter less than one in TRNG::chi_square_dist_pdf");
01133   if (x<=0.0)
01134     return 0.0;
01135   return std::exp((0.5*nu-1.0)*std::log(0.5*x)-
01136                   0.5*x-TRNG::ln_Gamma(0.5*nu))*0.5;
01137 }
01138 
01139 // ---------------------------------------------------------------------
01140 
01142 
01157 double TRNG::Student_t_dist_pdf(double x, double nu) {
01158   if (nu<=0.0)
01159     throw TRNG::error("parameter less than or equal zero in TRNG::Student_t_dist_pdf");
01160   const double pi=3.14159265358979324;
01161   double t=(nu+1.0)*0.5;
01162   return std::exp(TRNG::ln_Gamma(t)-TRNG::ln_Gamma(nu/2.0))/
01163     std::sqrt(pi*nu)*std::pow((1.0+x*x/nu), -t);
01164 }
01165 
01166 // ---------------------------------------------------------------------
01167 
01169 
01183 double TRNG::binomial_dist_pdf(long k, long n, double p=0.5) {
01184   if (p<=0.0 || p>1.0)
01185     throw TRNG::error("probability <=0.0 or >1.0 in TRNG::binomial_dist_pdf");
01186   if (n<1l)
01187     throw TRNG::error("less than one trail TRNG::binomial_dist_pdf");
01188   if (k<0 || k>n)
01189     return 0.0;
01190   return TRNG::binomial_coeff(n, k)*
01191     std::pow(p, static_cast<double>(k))*
01192     std::pow(1.0-p, static_cast<double>(n-k));
01193 }
01194 
01195 // ---------------------------------------------------------------------
01196 
01198 
01214 double TRNG::poisson_dist_pdf(long k, double mu) {
01215   if (k<0l)
01216     return 0.0;
01217   return std::pow(mu, static_cast<double>(k))*
01218     std::exp(-mu-TRNG::ln_factorial(k));
01219 }
01220 
01221 // ---------------------------------------------------------------------
01222 
01224 
01237 double TRNG::geometric_dist_pdf(long k, double q) {
01238   if (q<=0.0 ||q>1.0)
01239     throw TRNG::error("parameter out of range in TRNG::RNG::geometric_dist_pdf");
01240   if (k<=0l)
01241     return 0.0;
01242   if (k==1l)
01243     return q;
01244   return q*std::pow(1.0-q, static_cast<double>(k)-1.0);
01245 }

Generated on Wed Feb 19 02:00:03 2003 for Tina's Random Number Generators by doxygen1.2.15