00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023 #include "TRandom2.h"
00024 #include "TRandom3.h"
00025
00026
00027 ClassImp(TRandom2)
00028
00029
00030 TRandom2::TRandom2(UInt_t seed)
00031 {
00032
00033
00034
00035 SetName("Random2");
00036 SetTitle("Random number generator with period of about 10**26");
00037 SetSeed(seed);
00038 }
00039
00040
00041 TRandom2::~TRandom2()
00042 {
00043
00044
00045
00046 }
00047
00048
00049 Double_t TRandom2::Rndm(Int_t)
00050 {
00051
00052
00053
00054
00055 #define TAUSWORTHE(s,a,b,c,d) (((s &c) <<d) & 0xffffffffUL ) ^ ((((s <<a) & 0xffffffffUL )^s) >>b)
00056
00057 const double kScale = 2.3283064365386963e-10;
00058
00059 fSeed = TAUSWORTHE (fSeed, 13, 19, 4294967294UL, 12);
00060 fSeed1 = TAUSWORTHE (fSeed1, 2, 25, 4294967288UL, 4);
00061 fSeed2 = TAUSWORTHE (fSeed2, 3, 11, 4294967280UL, 17);
00062
00063 UInt_t iy = fSeed ^ fSeed1 ^ fSeed2;
00064 if (iy) return kScale*static_cast<Double_t>(iy);
00065 return Rndm();
00066 }
00067
00068
00069 void TRandom2::RndmArray(Int_t n, Float_t *array)
00070 {
00071
00072
00073 const double kScale = 2.3283064365386963e-10;
00074
00075 UInt_t iy;
00076
00077 for(Int_t i=0; i<n; i++) {
00078 fSeed = TAUSWORTHE (fSeed, 13, 19, 4294967294UL, 12);
00079 fSeed1 = TAUSWORTHE (fSeed1, 2, 25, 4294967288UL, 4);
00080 fSeed2 = TAUSWORTHE (fSeed2, 3, 11, 4294967280UL, 17);
00081
00082 iy = fSeed ^ fSeed1 ^ fSeed2;
00083 if (iy) array[i] = (Float_t)(kScale*static_cast<Double_t>(iy));
00084 else array[i] = Rndm();
00085 }
00086 }
00087
00088
00089 void TRandom2::RndmArray(Int_t n, Double_t *array)
00090 {
00091
00092
00093 const double kScale = 2.3283064365386963e-10;
00094
00095 UInt_t iy;
00096 for(Int_t i=0; i<n; i++) {
00097 fSeed = TAUSWORTHE (fSeed, 13, 19, 4294967294UL, 12);
00098 fSeed1 = TAUSWORTHE (fSeed1, 2, 25, 4294967288UL, 4);
00099 fSeed2 = TAUSWORTHE (fSeed2, 3, 11, 4294967280UL, 17);
00100
00101 iy = fSeed ^ fSeed1 ^ fSeed2;
00102 if (iy) array[i] = kScale*static_cast<Double_t>(iy);
00103 else array[i] = Rndm();
00104 }
00105 }
00106
00107
00108 void TRandom2::SetSeed(UInt_t seed)
00109 {
00110
00111
00112
00113
00114
00115
00116
00117
00118 #define LCG(n) ((69069 * n) & 0xffffffffUL) // linear congurential generator
00119
00120 if (seed > 0) {
00121 fSeed = LCG (seed);
00122 if (fSeed < 2) fSeed += 2UL;
00123 fSeed1 = LCG (fSeed);
00124 if (fSeed1 < 8) fSeed1 += 8UL;
00125 fSeed2 = LCG (fSeed1);
00126 if (fSeed2 < 16) fSeed2 += 16UL;
00127 } else {
00128
00129 TRandom3 r3(0);
00130 fSeed = static_cast<UInt_t> (4294967296.*r3.Rndm());
00131 fSeed1 = static_cast<UInt_t> (4294967296.*r3.Rndm());
00132 fSeed2 = static_cast<UInt_t> (4294967296.*r3.Rndm());
00133
00134 if (fSeed < 2) fSeed += 2UL;
00135 if (fSeed1 < 8) fSeed1 += 8UL;
00136 if (fSeed2 < 16) fSeed2 += 16UL;
00137 }
00138
00139
00140 for (int i = 0; i < 6; ++i)
00141 Rndm();
00142
00143 return;
00144 }
00145