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

trng.h

00001 // ---------------------------------------------------------------------
00002 // Time-stamp: <Tuesday, 10.12.2002, 13:14:42; edited by bauke>
00003 // 
00004 // Tina's random number generators TRNG
00005 //
00006 // Copyright (C) 2001, 2002 Heiko Bauke
00007 //
00008 // heiko.bauke@physik.uni-magdeburg.de
00009 //
00010 // TRNG is free software; you can redistribute it and/or
00011 // modify it under the terms of the GNU General Public License
00012 // as published by the Free Software Foundation. This program
00013 // is distributed WITHOUT ANY WARRANTY; without even the implied
00014 // warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
00015 // See the GNU General Public License for more details.
00016 //
00017 // ---------------------------------------------------------------------
00018 
00019 #if !defined TRNG_H
00020 
00021 #define TRNG_H
00022 
00023 #include <iostream>
00024 #include <cmath>
00025 #include <climits>
00026 #include <vector>
00027 #include <trnglib.h>
00028 
00029 #define TRNG_VERSION "2"
00030 #define TRNG_REVISION "1"
00031 
00032 
00034 
00037 namespace TRNG {
00038   
00040 
00044   enum RNG_type {
00045     RNG_t,          
00046     generic_MLCG_t, 
00047     ParkMiller_t,   
00048     LCG32_t,        
00049     LCG64_t,        
00050     MRG2_t,         
00051     MRG3_t,         
00052     MRG4_t,         
00053     CLCG2_t,        
00054     CLCG3_t,        
00055     CLCG4_t,        
00056     EINV_t,         
00057     EINVLCG64_t,    
00058     trng_gsl_t,     
00059     user1_t,        
00060     user2_t,        
00061     user3_t         
00062   };   
00063  
00065   struct vector2d_struct {
00067     double x1;
00069     double x2;
00070   };
00071 
00073 
00077   typedef struct vector2d_struct vector2d;
00078   
00080   struct vector3d_struct {
00082     double x1;
00084     double x2;
00086     double x3;
00087   };
00088 
00090 
00094   typedef struct vector3d_struct vector3d;
00095 
00097   struct vector4d_struct {
00099     double x1;
00101     double x2;
00103     double x3;
00105     double x4;
00106   };
00107 
00109 
00113   typedef struct vector4d_struct vector4d;
00114 
00115   // -------------------------------------------------------------------
00116 
00118 
00127   template<class RNG_type> class RNG {
00128   private:
00129     RNG_type & derived(void) {
00130       return static_cast<RNG_type &>(*this);
00131     }
00132     
00133   protected:
00135     long max_val;
00137     long max_val2;
00138   public:
00140 
00143     static const TRNG::RNG_type type=RNG_t;
00144     
00146 
00151     const char * name(void) {
00152       return derived().name(); 
00153     }
00154     
00156 
00160     void reset(void)  { 
00161       derived().reset(); 
00162     }
00163     
00165 
00173     void seed(long s=0l) { 
00174       derived().seed(s);
00175     }
00176     
00178 
00185     long rand(void) { 
00186       return derived().rand(); 
00187     }
00188 
00190 
00193     long max(void) {
00194       return max_val;
00195     };
00196     
00198 
00202     bool boolean(void) {
00203       return (rand()<max_val2) ? true : false;
00204     }
00205 
00207 
00211     bool boolean(const double p) {
00212       return (rand()<static_cast<long>(p*max())+1.0) ? true : false;
00213     }
00214     
00216 
00219     double uniform(void) {
00220       return uniformco();
00221     }
00222 
00224 
00229     double uniform(const double a, const double b) {
00230       return uniformco(a, b);
00231     }
00232 
00234 
00237     double uniformco(void) {
00238       return static_cast<double>(rand())/
00239         (static_cast<double>(max())+1.0);
00240     }
00241 
00243 
00248     double uniformco(const double a, const double b) {
00249       return a+(b-a)*static_cast<double>(rand())/
00250         (static_cast<double>(max())+1.0);
00251     }
00252 
00254 
00257     double uniformcc(void) {
00258       return static_cast<double>(rand())/
00259         static_cast<double>(max());
00260     }
00261     
00263 
00268     double uniformcc(const double a, const double b) {
00269       return a+(b-a)*static_cast<double>(rand())/
00270         static_cast<double>(max());
00271     }
00272 
00274 
00277     double uniformoc(void) {
00278       return (static_cast<double>(rand())+1.0)/
00279         (static_cast<double>(max())+1.0);
00280     }
00281     
00283 
00288     double uniformoc(const double a, const double b) {
00289       return a+(b-a)*(static_cast<double>(rand())+1.0)/
00290         (static_cast<double>(max())+1.0);
00291     }
00292     
00294 
00297     double uniformoo(void) {
00298       return (static_cast<double>(rand())+1.0)/
00299         (static_cast<double>(max())+2.0);
00300     }
00301 
00303 
00308     double uniformoo(const double a, const double b) {
00309       return a+(b-a)*(static_cast<double>(rand())+1.0)/
00310         (static_cast<double>(max())+2.0);
00311     }
00312     
00314 
00318     long uniforml(const long b) {
00319       return static_cast<long>(uniform(0.0, static_cast<double>(b)));
00320     }
00321     
00323 
00328     long uniforml(const long a, const long b) {
00329       return static_cast<long>(uniform(static_cast<double>(a),
00330                                        static_cast<double>(b)));
00331     }
00332     
00334 
00348     double normal_dist(const double sigma=1.0, const double mu=0.0) {
00349       // static double s, t1, t2;
00350       double s, t1, t2;
00351       // static bool calc=true;
00352       if (sigma<=0.0)
00353         throw error("negative or zero standard deviation in TRNG::RNG::normal_dist");
00354       // if (calc) {
00355       // calc=false;
00356       do {
00357         t1=uniformoo(-1.0, 1.0);
00358         t2=uniformoo(-1.0, 1.0);
00359         s=t1*t1+t2*t2;
00360       } while (s>=1.0);
00361       s=std::sqrt(-2.0*std::log(s)/s);
00362       return t1*s*sigma+mu;
00363       // } else {
00364       // calc=true;
00365       // return t2*s*sigma+mu;
00366       // }
00367     }
00368 
00370 
00389     double exp_dist(const double mu=1.0) {
00390       if (mu<=0.0)
00391         throw error("negative or zero parameter in TRNG::RNG::exp_dist");
00392       double t=uniformoc();
00393       return -mu*std::log(t);
00394     }
00395 
00397 
00411     double laplace_dist(const double a=1.0) {
00412       if (a<=0.0)
00413         throw error("parameter less or equal zero in TRNG::RNG::laplace_dist");
00414       double u;
00415       do {
00416         u=uniformoo(-1.0, 1.0);
00417       } while (u==0.0);
00418       return u<0.0 ? a*std::log(-u) : -a*std::log(u);
00419     }
00420 
00422 
00442     double tent_dist(const double a=1.0) {
00443       if (a<=0.0)
00444         throw error("parameter less or equal zero in TRNG::RNG::tent_dist");
00445       double y=uniformoo();
00446       if (y<0.5)
00447         return (-1.0+std::sqrt(2*y))*a;
00448       else
00449         return (1.0-std::sqrt(2.0-2.0*y))*a;
00450     }
00451     
00453 
00472     double Gamma_dist(const double a, const double b) {
00473       if (a<=0.0 || b<=0.0)
00474         throw error("parameter less or equal zero in TRNG::RNG::Gamma_dist");
00475       double a_int=std::floor(a);
00476       if (a_int==a) {
00477         // if a is an integer
00478         if (a<12) {
00479           double prod=1.0;
00480           for (int i=0; i<a; ++i)
00481             prod*=uniformoo();
00482           // Note: for 12 iterations we are safe against underflow, since
00483           // the smallest positive random number is O(2^-32). This means
00484           // the smallest possible product is 2^(-12*32) = 10^-116 which
00485           // is within the range of double precision.
00486           return -b*std::log(prod);
00487         } else {
00488           // Works only if a>1, and is most efficient if a is large
00489           // This algorithm, reported in Knuth, is attributed to Ahrens.  A
00490           // faster one, we are told, can be found in: J. H. Ahrens and
00491           // U. Dieter, Computing 12 (1974) 223-246.
00492           double sqa, x, y, v;
00493           sqa=std::sqrt(2.0*a-1.0);
00494           do {
00495             do {
00496               const double pi=3.14159265358979324;
00497               y=std::tan(pi*uniformco());
00498               x=sqa*y+a-1.0;
00499             }
00500             while (x<=0.0);
00501             v=uniformco();
00502           }
00503           while (v>(1.0+y*y)*std::exp((a-1.0)*log(x/(a-1))-sqa*y));
00504           return b*x;
00505         }
00506       }
00507       if (a_int==0.0) {
00508         // if a<1
00509         // This is exercise 16 from Knuth; see page 135, and the solution is
00510         // on page 551.
00511         const double e=2.71828182844;
00512         double p, q, x, u, v;
00513         p=e/(a+e);
00514         do {
00515           u=uniformco();
00516           v=uniformoo();
00517           if (u<p) {
00518             x=std::exp((1.0/a)*std::log(v));
00519             q=std::exp(-x);
00520           } else {
00521             x=1.0-std::log(v);
00522             q=std::exp((a-1.0)*std::log(x));
00523           }
00524         } while (uniformco()>=q);
00525         return b*x;
00526       }
00527       return Gamma_dist(a_int, b)+Gamma_dist(a-a_int, b); 
00528     }
00529 
00531 
00550     double Beta_dist(const double a, const double b) {
00551       if (a<=0.0 || b<=0.0)
00552         throw error("negative or zero parameter in TRNG::RNG::Beta_dist");
00553       double t1=Gamma_dist(a, 1.0);
00554       double t2=Gamma_dist(b, 1.0);
00555       return t1/(t1+t2);
00556     }
00557 
00559 
00575     double chi_square_dist(const double nu) {
00576       if (nu<1.0)
00577         throw error("parameter less than one in TRNG::RNG::chi_square_dist");
00578       return 2.0*Gamma_dist(0.5*nu, 1.0);
00579     }
00580 
00582 
00597     long binomial_dist(long n, double p=0.5) {
00598       if (p<=0.0 || p>1.0)
00599         throw error("probability <=0.0 or >1.0 in TRNG::RNG::binomial_dist");
00600       if (n<1l)
00601         throw error("less than one trail TRNG::RNG::binomial_dist");
00602       long i, a, b, k=0;
00603       while (n>12l) {  // This parameter is tunable
00604         double X;
00605         a=1l+(n/2l);
00606         b=1l+n-a;
00607         X=Beta_dist(static_cast<double>(a), static_cast<double>(b));
00608         if (X>=p) {
00609           n=a-1l;
00610           p/=X;
00611         } else {
00612           k+=a;
00613           n=b-1l;
00614           p=(p-X)/(1l-X);
00615         }
00616       }
00617       for (i=0l; i<n; i++)
00618         if (uniformco()<p)
00619           k++;
00620       return k;
00621     }
00622 
00624 
00637     long binomial_dist_tab(long n, double p=0.5) {
00638       if (p<=0.0 || p>1.0)
00639         throw error("probability <=0.0 or >1.0 in TRNG::RNG::binomial_dist_tab");
00640       if (n<1l)
00641         throw error("less than one trail TRNG::RNG::binomial_dist_tab");
00642       static long n_=0l;
00643       static double p_=0.0;
00644       static std::vector<double> prob;
00645       if (n!=n_ || p!=p_) {
00646         n_=n;
00647         p_=p;
00648         prob.resize(n_);
00649         // use
00650         // ( n )            (  n  )   n-k ( n )
00651         // (   ) = 1  and   (     ) = --- (   )
00652         // ( 0 )            ( k+1 )   k+1 ( k )
00653         double b=1.0;
00654         for (int k=0l; k<n_; ++k) {
00655           prob[k]=b*std::pow(p_, static_cast<double>(k))*
00656             std::pow(1.0-p_, static_cast<double>(n_-k));
00657           if (k>0l)
00658             prob[k]+=prob[k-1l];
00659           b=std::floor(b*static_cast<double>(n_-k)/
00660                        static_cast<double>(k+1l)+0.5); 
00661         }
00662       }
00663       return find_interval(prob, uniformco());
00664     }
00665 
00667 
00683     double Student_t_dist(const double nu) {
00684       if (nu<=0.0)
00685         throw error("parameter less than or equal zero in TRNG::RNG::Student_t_dist");
00686       if (nu<=2.0) 
00687         return normal_dist()/std::sqrt(chi_square_dist(nu)/nu);
00688       else {
00689         double y1, y2, z;
00690         do {
00691           y1=normal_dist();
00692           y2=exp_dist(1.0/(0.5*nu-1.0));
00693           z=y1*y1/(nu-2.0);
00694         } while (1.0-z<0.0 || std::exp(-y2-z)>(1.0-z));
00695         // Note that there is a typo in Knuth's formula, the line below
00696         // is taken from the original paper of Marsaglia, Mathematics of
00697         // Computation, 34 (1980), p 234-256 
00698         return y1/std::sqrt((1.0-2.0/nu)*(1.0-z));
00699       }
00700     }
00701 
00703 
00715     long poisson_dist(double mu=1.0) {
00716       if (mu<=0)
00717         throw error("parameter less than 0 in TRNG::RNG::poisson_dist");
00718       double emu;
00719       double prod=1.0;
00720       long k=0;  
00721       while (mu>10.0) {
00722         long m=static_cast<long>(mu*(7.0/8.0));
00723         double X=Gamma_dist(static_cast<long>(m), 1.0);
00724         if (X>=mu) {
00725           return k+binomial_dist(m-1l, mu/X);
00726         } else {
00727           k+=m;
00728           mu-=X; 
00729         }
00730       }
00731       // This following method works well when mu is small
00732       emu=std::exp(-mu);
00733       do {
00734         prod*=uniformco();
00735         k++;
00736       } while (prod>emu);
00737       return k-1l;
00738     }
00739 
00741 
00753     long geometric_dist(double q) {
00754       if (q<=0.0 ||q>1.0)
00755         throw error("parameter out of range in TRNG::RNG::geometric_dist");
00756       if (q==1.0)
00757         return 1l;
00758       return static_cast<long>
00759         (std::log(uniformoo())/std::log(1.0-p))+1l;
00760     }
00761 
00763 
00789     template<class t_function>
00790       double rejection(t_function p, double a1, double a2, double p_max)  {
00791       double t1, t2;
00792       do {
00793         t1=uniformco(a1, a2);
00794         t2=uniformco(0.0, p_max);
00795       } while (p(t1)<t2);
00796       return t1;
00797     }
00798 
00800 
00812     long discrete_dist(const std::vector<double> p) {
00813       if (p.empty())
00814         throw error("empty vector in TRNG::discrete_dist");
00815       double x=uniformco();
00816       if (x<p[0])
00817         return 0l;
00818       long i1=0l;
00819       long i2=p.size()-1l;
00820       while (i2-i1>1l) {
00821         long i3=(i2+i1)/2;
00822         if (x<=p[i3])
00823           i2=i3;
00824         else
00825           i1=i3;
00826       }
00827       return i2;
00828     }
00829 
00831 
00839     vector2d spherical2d(void) {
00840       double r2, t1, t2, t12, t22;
00841       vector2d vec;
00842       do {
00843         t1=2.0*uniformco()-1.0;
00844         t2=uniformco();
00845         t12=t1*t1;
00846         t22=t2*t2;
00847         r2=t12+t22;
00848       } while (r2>1.0);
00849       vec.x1=2.0*t1*t2/r2;
00850       vec.x2=(t12-t22)/r2;
00851       return vec;
00852     }
00853 
00855 
00863     vector3d spherical3d(void)  {
00864       double q, r2, t1, t2;
00865       vector3d vec;
00866       do {
00867         t1=2.0*uniformco()-1.0;
00868         t2=2.0*uniformco()-1.0;
00869         r2=t1*t1+t2*t2;
00870       } while (r2>1.0);
00871       q=2.0*sqrt(1.0-r2);
00872       vec.x1=t1*q;
00873       vec.x2=t2*q;
00874       vec.x3=1-2*r2;
00875       return vec;
00876     }
00877 
00879 
00887     vector4d spherical4d(void)  {
00888       double q, r21, r22, t1, t2, t3, t4;
00889       vector4d vec;
00890       do {
00891         t1=2.0*uniformco()-1.0;
00892         t2=2.0*uniformco()-1.0;
00893         r21=t1*t1+t2*t2;
00894       } while (r21>1.0);
00895       do {
00896         t3=2.0*uniformco()-1.0;
00897         t4=2.0*uniformco()-1.0;
00898         r22=t3*t3+t4*t4;
00899       } while (r22>1.0);
00900       q=sqrt((1.0-r21)/r22);
00901       vec.x1=t1;
00902       vec.x2=t2;
00903       vec.x3=t3*q;
00904       vec.x4=t4*q;
00905       return(vec);
00906     }
00907 
00909 
00918     void split(long s, long n) {
00919       derived().split(s, n);
00920     }
00921 
00923 
00929     void jump(long long s) {
00930       if (s<0ll)
00931         throw error("invalid argument for TRNG::RNG::jump");
00932       long i=0l;
00933       while (s>0ll) {
00934         if (s%2==1)
00935           jump2(i);
00936         ++i;
00937         s>>=1;
00938       }
00939     }
00940 
00942 
00948     void jump(long long s, long n) {
00949       if (s<0ll || n<0l)
00950         throw error("invalid argument for TRNG::RNG::jump");
00951       while (n>0l) {
00952         jump(s);
00953         --n;
00954       }
00955     }
00956 
00958 
00964     void jump2(long s) {
00965       derived().jump2(s);
00966     }
00967 
00969 
00976     void jump2(long s, long n) {
00977       if (s<0l || n<0l)
00978         throw error("invalid argument for TRNG::RNG::jump2");
00979       while (n>0l) {
00980         if (n%2l==1l)
00981           jump2(s);
00982         s+=1l;
00983         n>>=1;
00984       }
00985     }
00986 
00988 
00993     void save_status(std::vector<long> &s) {
00994       derived().save_status(s);
00995     }
00996 
00998 
01003     void load_status(const std::vector<long> &s) {
01004       derived().load_status(s);
01005     }
01006 
01007 //     //! generic copy constructor
01008 //     /*!
01009 //       This is the random number generator copy constructor. Useful if
01010 //       some generators with the same sequence of random numbers are needed.
01011 //      */
01012 //     RNG<RNG_type> & RNG<RNG_type>::operator=(RNG<RNG_type> &other) {
01013 //       if (this!=&other) {
01014 //      std::vector<long> s;
01015 //      other.save_status(s);
01016 //      load_status(s);
01017 //       }
01018 //       return *this;
01019 //     }
01020     
01021   };
01022 
01023   // --------------------------------------------------------------------
01024 
01025   class generic_MLCG : public RNG<generic_MLCG> {
01026     // generic class for multiplicative linear congruential 
01027     // random number generators
01028   private:
01029     long a, a_save;
01030     long modulus, modulus_save;
01031     long r;
01032     void backward(void);
01033   public:
01034     static const RNG_type type=generic_MLCG_t;
01035     const char * name(void);    
01036     void reset(void);
01037     void seed(long s=0l);
01038     long rand(void);
01039     void split(long s, long n);    
01040     void jump2(long s);
01041     void save_status(std::vector<long> &s);
01042     void load_status(const std::vector<long> &s);
01043     generic_MLCG & generic_MLCG::operator=(RNG<generic_MLCG> &other);
01044     generic_MLCG(long a_=16807l, long modulus_=2147483647l, long seed_=0l);
01045   };
01046 
01047   // -------------------------------------------------------------------
01048 
01050 
01065   class ParkMiller : public RNG<ParkMiller> {
01066   private:
01067     generic_MLCG R;
01068   public:
01069     static const RNG_type type=ParkMiller_t;
01070     const char * name(void);    
01071     void reset(void);
01072     void seed(long s=0l);
01073     long rand(void);    
01074     void split(long s, long n);
01075     void jump2(long s);
01076     void save_status(std::vector<long> &s);
01077     void load_status(const std::vector<long> &s);
01078     ParkMiller & ParkMiller::operator=(RNG<ParkMiller> &other); 
01080 
01089     ParkMiller(long a_=16807l, long modulus_=2147483647l, long seed_=0l);
01090   };
01091   
01092   // -------------------------------------------------------------------
01093 
01095 
01109   class LCG32 : public RNG<LCG32> {
01110   private:
01111     unsigned long a, a_save;
01112     unsigned long b, b_save;
01113     unsigned long r;
01114     void backward(void);
01115   public:
01116     static const RNG_type type=LCG32_t;
01117     const char * name(void);
01118     void reset(void);
01119     void seed(long s=0l);    
01120     long rand(void);    
01121     void split(long s, long n);    
01122     void jump2(long s);    
01123     void save_status(std::vector<long> &s);    
01124     void load_status(const std::vector<long> &s);
01125     LCG32 & LCG32::operator=(RNG<LCG32> &other);
01127 
01136     LCG32(unsigned long a_=69069ul, unsigned long b_=1ul, long seed_=0l);
01137   };
01138 
01139   // -------------------------------------------------------------------
01140 
01142 
01156   class LCG64 : public RNG<LCG64> {
01157   private:
01158     unsigned long long a, a_save;
01159     unsigned long long b, b_save;
01160     unsigned long long r;
01161     void backward(void);
01162   public:
01163     static const RNG_type type=LCG64_t;
01164     const char * name(void);
01165     void reset(void);
01166     void seed(long s=0l);
01167     long rand(void);
01168     void split(long s, long n);
01169     void jump2(long s);    
01170     void save_status(std::vector<long> &s);    
01171     void load_status(const std::vector<long> &s);
01172     LCG64 & LCG64::operator=(RNG<LCG64> &other);
01174 
01183     LCG64(unsigned long long a_=18145460002477866997ull, 
01184           unsigned long long b_=1ul, 
01185           long seed_=0ul);    
01186   };
01187   
01188   // -------------------------------------------------------------------
01189 
01191 
01199   class MRG2 : public RNG<MRG2> {
01200   private:
01201     long a0, a1, a0_save, a1_save;
01202     long r0, r1;
01203     long modulus, modulus_save;
01204     void backward(void);
01205   public:
01206     static const RNG_type type=MRG2_t;
01207     const char * name(void);    
01208     void reset(void);
01209     void seed(long s=0l);
01210     long rand(void);    
01211     void split(long s, long n);
01212     void jump2(long s);
01213     void save_status(std::vector<long> &s);
01214     void load_status(const std::vector<long> &s);    
01215     MRG2 & MRG2::operator=(RNG<MRG2> &other);
01217 
01229     MRG2(long a0_=1498809829l, long a1_=1160990996l, 
01230          long modulus_=2147483647l, long seed_=0l);
01231   };
01232 
01233   // -------------------------------------------------------------------
01234 
01236 
01244   class MRG3 : public RNG<MRG3> {
01245   private:
01246     long a0, a1, a2, a0_save, a1_save, a2_save;
01247     long r0, r1, r2;
01248     long modulus, modulus_save;
01249     void backward(void);
01250   public:
01251     static const RNG_type type=MRG3_t;
01252     const char * name(void);
01253     void reset(void);    
01254     void seed(long s=0l);
01255     long rand(void);   
01256     void split(long s, long n);
01257     void jump2(long s);    
01258     void save_status(std::vector<long> &s);
01259     void load_status(const std::vector<long> &s);
01260     MRG3 & MRG3::operator=(RNG<MRG3> &other);
01262 
01275     MRG3(long a0_=2021422057l, long a1_=1826992351l, 
01276          long a2_=1977753457l, long modulus_=2147483647l, 
01277          long seed_=0l);
01278   };
01279     
01280   // -------------------------------------------------------------------
01281 
01283 
01292   class MRG4 : public RNG<MRG4> {
01293   private:
01294     long a0, a1, a2, a3, a0_save, a1_save, a2_save, a3_save;
01295     long r0, r1, r2, r3;
01296     long modulus, modulus_save;
01297     void backward(void);
01298   public:
01299     static const RNG_type type=MRG4_t;
01300     const char * name(void);
01301     void reset(void);
01302     void seed(long s=0l);
01303     long rand(void);
01304     void split(long s, long n);
01305     void jump2(long s);
01306     void save_status(std::vector<long> &s);
01307     void load_status(const std::vector<long> &s);
01308     MRG4 & MRG4::operator=(RNG<MRG4> &other);
01310 
01325     MRG4(long a0_=2001982722l, long a1_=1412284257l, 
01326          long a2_=1155380217l, long a3_=1668339922l,
01327          long modulus_=2147483647l, long seed_=0l);
01328   };
01329 
01330   // -------------------------------------------------------------------
01331 
01333 
01346   class CLCG2 : public RNG<CLCG2> {
01347   private:
01348     generic_MLCG R1;
01349     generic_MLCG R2;
01350     long modulus;
01351   public:
01352     static const RNG_type type=CLCG2_t;
01353     const char * name(void);
01354     void reset(void);
01355     void seed(long s=0l);
01356     long rand(void);    
01357     void split(long s, long n);
01358     void jump2(long s);
01359     void save_status(std::vector<long> &s);    
01360     void load_status(const std::vector<long> &s);    
01361     CLCG2 & CLCG2::operator=(RNG<CLCG2> &other);
01363 
01383     CLCG2(long a1= 376555083l, long m1=2147482951l,
01384           long a2=1028879659l, long m2=2147482949l, 
01385           long seed_=0l);    
01386   };
01387 
01388   // -------------------------------------------------------------------
01389 
01391 
01405   class CLCG3 : public RNG<CLCG3> {
01406   private:
01407     generic_MLCG R1;
01408     generic_MLCG R2;
01409     generic_MLCG R3;
01410     long modulus;
01411   public:
01412     static const RNG_type type=CLCG3_t;
01413     const char * name(void);
01414     void reset(void);
01415     void seed(long s=0l);
01416     long rand(void);
01417     void split(long s, long n);
01418     void jump2(long s);
01419     void save_status(std::vector<long> &s);
01420     void load_status(const std::vector<long> &s);
01421     CLCG3 & CLCG3::operator=(RNG<CLCG3> &other);
01423 
01447     CLCG3(long a1= 376555083l, long m1=2147482951l,
01448           long a2=1028879659l, long m2=2147482949l,
01449           long a3= 225802979l, long m3=2147482943l,
01450           long seed_=0l);    
01451   };
01452 
01453   // -------------------------------------------------------------------
01454 
01456 
01471   class CLCG4 : public RNG<CLCG4> {
01472   private:
01473     generic_MLCG R1;
01474     generic_MLCG R2;
01475     generic_MLCG R3;
01476     generic_MLCG R4;
01477     long modulus;
01478   public:
01479     static const RNG_type type=CLCG4_t;
01480     const char * name(void);
01481     void reset(void);
01482     void seed(long s=0l);
01483     long rand(void);
01484     void split(long s, long n);
01485     void jump2(long s);
01486     void save_status(std::vector<long> &s);  
01487     void load_status(const std::vector<long> &s);
01488     CLCG4 & CLCG4::operator=(RNG<CLCG4> &other);
01490 
01518     CLCG4(long a1= 376555083l, long m1=2147482951l,
01519           long a2=1028879659l, long m2=2147482949l,
01520           long a3= 225802979l, long m3=2147482943l,
01521           long a4=2028073966l, long m4=2147482859l,
01522           long seed_=0l);    
01523   };
01524 
01525   // -------------------------------------------------------------------
01526 
01528 
01541   class EINV : public RNG<EINV> {
01542   private:
01543     long a, i, di, modulus, a_save;
01544     long modulo_invers(const long a);
01545   public:
01546     static const RNG_type type=EINV_t;
01547     const char * name(void);
01548     void reset(void);
01549     void seed(long s=0l);
01550     long rand(void);
01551     void split(long s, long n);
01552     void jump2(long s);
01553     void save_status(std::vector<long> &s);
01554     void load_status(const std::vector<long> &s);
01555     EINV & EINV::operator=(RNG<EINV> &other);
01557 
01562     EINV(long a_=1073741831l, long seed_=0l);
01563     
01564   };
01565 
01566   // -------------------------------------------------------------------
01567 
01569 
01581   class EINVLCG64 : public RNG<EINVLCG64> {
01582   private:
01583     LCG64 R1;
01584     EINV R2;
01585   public:
01586     static const RNG_type type=EINVLCG64_t;
01587     const char * name(void);
01588     void reset(void);
01589     void seed(long s=0l);    
01590     long rand(void);
01591     void split(long s, long n);
01592     void jump2(long s);
01593     void save_status(std::vector<long> &s);
01594     void load_status(const std::vector<long> &s);    
01595     EINVLCG64 & EINVLCG64::operator=(RNG<EINVLCG64> &other);
01597 
01600     EINVLCG64(long seed_=0l);    
01601   };
01602 
01603 }
01604 
01605 
01606 #endif

Generated on Tue Dec 10 13:31:37 2002 for Tina's Random Number Generators by doxygen1.2.15