00001 // @(#)root/mathcore:$Id: TRandom.cxx 33362 2010-05-04 08:00:22Z moneta $ 00002 // Author: Rene Brun, Lorenzo Moneta 15/12/95 00003 00004 /************************************************************************* 00005 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. * 00006 * All rights reserved. * 00007 * * 00008 * For the licensing terms see $ROOTSYS/LICENSE. * 00009 * For the list of contributors see $ROOTSYS/README/CREDITS. * 00010 *************************************************************************/ 00011 00012 ////////////////////////////////////////////////////////////////////////// 00013 // 00014 // TRandom 00015 // 00016 // basic Random number generator class (periodicity = 10**9). 00017 // Note that this is a very simple generator (linear congruential) 00018 // which is known to have defects (the lower random bits are correlated) 00019 // and therefore should NOT be used in any statistical study. 00020 // One should use instead TRandom1, TRandom2 or TRandom3. 00021 // TRandom3, is based on the "Mersenne Twister generator", and is the recommended one, 00022 // since it has good random proprieties (period of about 10**6000 ) and it is fast. 00023 // TRandom1, based on the RANLUX algorithm, has mathematically proven random proprieties 00024 // and a period of about 10**171. It is however slower than the others. 00025 // TRandom2, is based on the Tausworthe generator of L'Ecuyer, and it has the advantage 00026 // of being fast and using only 3 words (of 32 bits) for the state. The period is 10**26. 00027 // 00028 // The following table shows some timings (in nanoseconds/call) 00029 // for the random numbers obtained using an Intel Pentium 3.0 GHz running Linux 00030 // and using the gcc 3.2.3 compiler 00031 // 00032 // TRandom 34 ns/call (BAD Generator) 00033 // TRandom1 242 ns/call 00034 // TRandom2 37 ns/call 00035 // TRandom3 45 ns/call 00036 // 00037 // 00038 // The following basic Random distributions are provided: 00039 // =================================================== 00040 // -Exp(tau) 00041 // -Integer(imax) 00042 // -Gaus(mean,sigma) 00043 // -Rndm() 00044 // -Uniform(x1) 00045 // -Landau(mpv,sigma) 00046 // -Poisson(mean) 00047 // -Binomial(ntot,prob) 00048 // 00049 // Random numbers distributed according to 1-d, 2-d or 3-d distributions 00050 // ===================================================================== 00051 // contained in TF1, TF2 or TF3 objects. 00052 // For example, to get a random number distributed following abs(sin(x)/x)*sqrt(x) 00053 // you can do : 00054 // TF1 *f1 = new TF1("f1","abs(sin(x)/x)*sqrt(x)",0,10); 00055 // double r = f1->GetRandom(); 00056 // or you can use the UNURAN package. You need in this case to initialize UNURAN 00057 // to the function you would like to generate. 00058 // TUnuran u; 00059 // u.Init(TUnuranDistrCont(f1)); 00060 // double r = u.Sample(); 00061 // 00062 // The techniques of using directly a TF1,2 or 3 function is powerful and 00063 // can be used to generate numbers in the defined range of the function. 00064 // Getting a number from a TF1,2,3 function is also quite fast. 00065 // UNURAN is a powerful and flexible tool which containes various methods for 00066 // generate random numbers for continuous distributions of one and multi-dimension. 00067 // It requires some set-up (initialization) phase and can be very fast when the distribution 00068 // parameters are not changed for every call. 00069 // 00070 // The following table shows some timings (in nanosecond/call) 00071 // for basic functions, TF1 functions and using UNURAN obtained running 00072 // the tutorial math/testrandom.C 00073 // Numbers have been obtained on an Intel Xeon Quad-core Harpertown (E5410) 2.33 GHz running 00074 // Linux SLC4 64 bit and compiled with gcc 3.4 00075 // 00076 // Distribution nanoseconds/call 00077 // TRandom TRandom1 TRandom2 TRandom3 00078 // Rndm.............. 5.000 105.000 7.000 10.000 00079 // RndmArray......... 4.000 104.000 6.000 9.000 00080 // Gaus.............. 36.000 180.000 40.000 48.000 00081 // Rannor............ 118.000 220.000 120.000 124.000 00082 // Landau............ 22.000 123.000 26.000 31.000 00083 // Exponential....... 93.000 198.000 98.000 104.000 00084 // Binomial(5,0.5)... 30.000 548.000 46.000 65.000 00085 // Binomial(15,0.5).. 75.000 1615.000 125.000 178.000 00086 // Poisson(3)........ 96.000 494.000 109.000 125.000 00087 // Poisson(10)....... 138.000 1236.000 165.000 203.000 00088 // Poisson(70)....... 818.000 1195.000 835.000 844.000 00089 // Poisson(100)...... 837.000 1218.000 849.000 864.000 00090 // GausTF1........... 83.000 180.000 87.000 88.000 00091 // LandauTF1......... 80.000 180.000 83.000 86.000 00092 // GausUNURAN........ 40.000 139.000 41.000 44.000 00093 // PoissonUNURAN(10). 85.000 271.000 92.000 102.000 00094 // PoissonUNURAN(100) 62.000 256.000 69.000 78.000 00095 // 00096 // Note that the time to generate a number from an arbitrary TF1 function 00097 // using TF1::GetRandom or using TUnuran is independent of the complexity of the function. 00098 // 00099 // TH1::FillRandom(TH1 *) or TH1::FillRandom(const char *tf1name) 00100 // ============================================================== 00101 // can be used to fill an histogram (1-d, 2-d, 3-d from an existing histogram 00102 // or from an existing function. 00103 // 00104 // Note this interesting feature when working with objects 00105 // ======================================================= 00106 // You can use several TRandom objects, each with their "independent" 00107 // random sequence. For example, one can imagine 00108 // TRandom *eventGenerator = new TRandom(); 00109 // TRandom *tracking = new TRandom(); 00110 // eventGenerator can be used to generate the event kinematics. 00111 // tracking can be used to track the generated particles with random numbers 00112 // independent from eventGenerator. 00113 // This very interesting feature gives the possibility to work with simple 00114 // and very fast random number generators without worrying about 00115 // random number periodicity as it was the case with Fortran. 00116 // One can use TRandom::SetSeed to modify the seed of one generator. 00117 // 00118 // a TRandom object may be written to a Root file 00119 // ============================================== 00120 // -as part of another object 00121 // -or with its own key (example gRandom->Write("Random"); 00122 // 00123 ////////////////////////////////////////////////////////////////////////// 00124 00125 #include "TROOT.h" 00126 #include "TMath.h" 00127 #include "TRandom.h" 00128 #include "TRandom3.h" 00129 #include "TSystem.h" 00130 #include "TDirectory.h" 00131 #include "Math/QuantFuncMathCore.h" 00132 #include <time.h> 00133 00134 ClassImp(TRandom) 00135 00136 //______________________________________________________________________________ 00137 TRandom::TRandom(UInt_t seed): TNamed("Random","Default Random number generator") 00138 { 00139 //*-*-*-*-*-*-*-*-*-*-*default constructor*-*-*-*-*-*-*-*-*-*-*-*-*-*-* 00140 //*-* =================== 00141 00142 SetSeed(seed); 00143 00144 } 00145 00146 //______________________________________________________________________________ 00147 TRandom::~TRandom() 00148 { 00149 //*-*-*-*-*-*-*-*-*-*-*default destructor*-*-*-*-*-*-*-*-*-*-*-*-*-*-* 00150 //*-* ================== 00151 00152 if (gRandom == this) gRandom = 0; 00153 } 00154 00155 //______________________________________________________________________________ 00156 Int_t TRandom::Binomial(Int_t ntot, Double_t prob) 00157 { 00158 // Generates a random integer N according to the binomial law 00159 // Coded from Los Alamos report LA-5061-MS 00160 // 00161 // N is binomially distributed between 0 and ntot inclusive 00162 // with mean prob*ntot. 00163 // prob is between 0 and 1. 00164 // 00165 // Note: This function should not be used when ntot is large (say >100). 00166 // The normal approximation is then recommended instead 00167 // (with mean =*ntot+0.5 and standard deviation sqrt(ntot*prob*(1-prob)). 00168 00169 if (prob < 0 || prob > 1) return 0; 00170 Int_t n = 0; 00171 for (Int_t i=0;i<ntot;i++) { 00172 if (Rndm() > prob) continue; 00173 n++; 00174 } 00175 return n; 00176 } 00177 00178 //______________________________________________________________________________ 00179 Double_t TRandom::BreitWigner(Double_t mean, Double_t gamma) 00180 { 00181 // Return a number distributed following a BreitWigner function with mean and gamma 00182 00183 Double_t rval, displ; 00184 rval = 2*Rndm() - 1; 00185 displ = 0.5*gamma*TMath::Tan(rval*TMath::PiOver2()); 00186 00187 return (mean+displ); 00188 } 00189 00190 //______________________________________________________________________________ 00191 void TRandom::Circle(Double_t &x, Double_t &y, Double_t r) 00192 { 00193 // generates random vectors, uniformly distributed over a circle of given radius. 00194 // Input : r = circle radius 00195 // Output: x,y a random 2-d vector of length r 00196 00197 Double_t phi = Uniform(0,TMath::TwoPi()); 00198 x = r*TMath::Cos(phi); 00199 y = r*TMath::Sin(phi); 00200 } 00201 00202 //______________________________________________________________________________ 00203 Double_t TRandom::Exp(Double_t tau) 00204 { 00205 // returns an exponential deviate. 00206 // 00207 // exp( -t/tau ) 00208 00209 Double_t x = Rndm(); // uniform on ] 0, 1 ] 00210 Double_t t = -tau * TMath::Log( x ); // convert to exponential distribution 00211 return t; 00212 } 00213 00214 //______________________________________________________________________________ 00215 Double_t TRandom::Gaus(Double_t mean, Double_t sigma) 00216 { 00217 // 00218 // samples a random number from the standard Normal (Gaussian) Distribution 00219 // with the given mean and sigma. 00220 // Uses the Acceptance-complement ratio from W. Hoermann and G. Derflinger 00221 // This is one of the fastest existing method for generating normal random variables. 00222 // It is a factor 2/3 faster than the polar (Box-Muller) method used in the previous 00223 // version of TRandom::Gaus. The speed is comparable to the Ziggurat method (from Marsaglia) 00224 // implemented for example in GSL and available in the MathMore library. 00225 // 00226 // 00227 // REFERENCE: - W. Hoermann and G. Derflinger (1990): 00228 // The ACR Method for generating normal random variables, 00229 // OR Spektrum 12 (1990), 181-185. 00230 // 00231 // Implementation taken from 00232 // UNURAN (c) 2000 W. Hoermann & J. Leydold, Institut f. Statistik, WU Wien 00233 /////////////////////////////////////////////////////////////////////////////// 00234 00235 00236 00237 const Double_t kC1 = 1.448242853; 00238 const Double_t kC2 = 3.307147487; 00239 const Double_t kC3 = 1.46754004; 00240 const Double_t kD1 = 1.036467755; 00241 const Double_t kD2 = 5.295844968; 00242 const Double_t kD3 = 3.631288474; 00243 const Double_t kHm = 0.483941449; 00244 const Double_t kZm = 0.107981933; 00245 const Double_t kHp = 4.132731354; 00246 const Double_t kZp = 18.52161694; 00247 const Double_t kPhln = 0.4515827053; 00248 const Double_t kHm1 = 0.516058551; 00249 const Double_t kHp1 = 3.132731354; 00250 const Double_t kHzm = 0.375959516; 00251 const Double_t kHzmp = 0.591923442; 00252 /*zhm 0.967882898*/ 00253 00254 const Double_t kAs = 0.8853395638; 00255 const Double_t kBs = 0.2452635696; 00256 const Double_t kCs = 0.2770276848; 00257 const Double_t kB = 0.5029324303; 00258 const Double_t kX0 = 0.4571828819; 00259 const Double_t kYm = 0.187308492 ; 00260 const Double_t kS = 0.7270572718 ; 00261 const Double_t kT = 0.03895759111; 00262 00263 Double_t result; 00264 Double_t rn,x,y,z; 00265 00266 00267 do { 00268 y = Rndm(); 00269 00270 if (y>kHm1) { 00271 result = kHp*y-kHp1; break; } 00272 00273 else if (y<kZm) { 00274 rn = kZp*y-1; 00275 result = (rn>0) ? (1+rn) : (-1+rn); 00276 break; 00277 } 00278 00279 else if (y<kHm) { 00280 rn = Rndm(); 00281 rn = rn-1+rn; 00282 z = (rn>0) ? 2-rn : -2-rn; 00283 if ((kC1-y)*(kC3+TMath::Abs(z))<kC2) { 00284 result = z; break; } 00285 else { 00286 x = rn*rn; 00287 if ((y+kD1)*(kD3+x)<kD2) { 00288 result = rn; break; } 00289 else if (kHzmp-y<exp(-(z*z+kPhln)/2)) { 00290 result = z; break; } 00291 else if (y+kHzm<exp(-(x+kPhln)/2)) { 00292 result = rn; break; } 00293 } 00294 } 00295 00296 while (1) { 00297 x = Rndm(); 00298 y = kYm * Rndm(); 00299 z = kX0 - kS*x - y; 00300 if (z>0) 00301 rn = 2+y/x; 00302 else { 00303 x = 1-x; 00304 y = kYm-y; 00305 rn = -(2+y/x); 00306 } 00307 if ((y-kAs+x)*(kCs+x)+kBs<0) { 00308 result = rn; break; } 00309 else if (y<x+kT) 00310 if (rn*rn<4*(kB-log(x))) { 00311 result = rn; break; } 00312 } 00313 } while(0); 00314 00315 00316 return mean + sigma * result; 00317 00318 } 00319 00320 00321 00322 //______________________________________________________________________________ 00323 UInt_t TRandom::Integer(UInt_t imax) 00324 { 00325 // returns a random integer on [ 0, imax-1 ]. 00326 00327 UInt_t ui; 00328 ui = (UInt_t)(imax*Rndm()); 00329 return ui; 00330 } 00331 00332 //______________________________________________________________________________ 00333 Double_t TRandom::Landau(Double_t mpv, Double_t sigma) 00334 { 00335 // Generate a random number following a Landau distribution 00336 // with mpv(most probable value) and sigma 00337 // Use function landau_quantile(x,sigma) which provides 00338 // the inverse of the landau cumulative distribution 00339 // landau_quantile has been converted from CERNLIB ranlan(G110) 00340 00341 if (sigma <= 0) return 0; 00342 Double_t x = Rndm(); 00343 Double_t res = mpv + ROOT::Math::landau_quantile(x, sigma); 00344 return res; 00345 } 00346 00347 //______________________________________________________________________________ 00348 Int_t TRandom::Poisson(Double_t mean) 00349 { 00350 // Generates a random integer N according to a Poisson law. 00351 // Prob(N) = exp(-mean)*mean^N/Factorial(N) 00352 // 00353 // Use a different procedure according to the mean value. 00354 // The algorithm is the same used by CLHEP 00355 // For lower value (mean < 25) use the rejection method based on 00356 // the exponential 00357 // For higher values use a rejection method comparing with a Lorentzian 00358 // distribution, as suggested by several authors 00359 // This routine since is returning 32 bits integer will not work for values larger than 2*10**9 00360 // One should then use the Trandom::PoissonD for such large values 00361 // 00362 Int_t n; 00363 if (mean <= 0) return 0; 00364 if (mean < 25) { 00365 Double_t expmean = TMath::Exp(-mean); 00366 Double_t pir = 1; 00367 n = -1; 00368 while(1) { 00369 n++; 00370 pir *= Rndm(); 00371 if (pir <= expmean) break; 00372 } 00373 return n; 00374 } 00375 // for large value we use inversion method 00376 else if (mean < 1E9) { 00377 Double_t em, t, y; 00378 Double_t sq, alxm, g; 00379 Double_t pi = TMath::Pi(); 00380 00381 sq = TMath::Sqrt(2.0*mean); 00382 alxm = TMath::Log(mean); 00383 g = mean*alxm - TMath::LnGamma(mean + 1.0); 00384 00385 do { 00386 do { 00387 y = TMath::Tan(pi*Rndm()); 00388 em = sq*y + mean; 00389 } while( em < 0.0 ); 00390 00391 em = TMath::Floor(em); 00392 t = 0.9*(1.0 + y*y)* TMath::Exp(em*alxm - TMath::LnGamma(em + 1.0) - g); 00393 } while( Rndm() > t ); 00394 00395 return static_cast<Int_t> (em); 00396 00397 } 00398 else { 00399 // use Gaussian approximation vor very large values 00400 n = Int_t(Gaus(0,1)*TMath::Sqrt(mean) + mean +0.5); 00401 return n; 00402 } 00403 } 00404 00405 //______________________________________________________________________________ 00406 Double_t TRandom::PoissonD(Double_t mean) 00407 { 00408 // Generates a random number according to a Poisson law. 00409 // Prob(N) = exp(-mean)*mean^N/Factorial(N) 00410 // 00411 // This function is a variant of TRandom::Poisson returning a double 00412 // instead of an integer. 00413 // 00414 Int_t n; 00415 if (mean <= 0) return 0; 00416 if (mean < 25) { 00417 Double_t expmean = TMath::Exp(-mean); 00418 Double_t pir = 1; 00419 n = -1; 00420 while(1) { 00421 n++; 00422 pir *= Rndm(); 00423 if (pir <= expmean) break; 00424 } 00425 return static_cast<Double_t>(n); 00426 } 00427 // for large value we use inversion method 00428 else if (mean < 1E9) { 00429 Double_t em, t, y; 00430 Double_t sq, alxm, g; 00431 Double_t pi = TMath::Pi(); 00432 00433 sq = TMath::Sqrt(2.0*mean); 00434 alxm = TMath::Log(mean); 00435 g = mean*alxm - TMath::LnGamma(mean + 1.0); 00436 00437 do { 00438 do { 00439 y = TMath::Tan(pi*Rndm()); 00440 em = sq*y + mean; 00441 } while( em < 0.0 ); 00442 00443 em = TMath::Floor(em); 00444 t = 0.9*(1.0 + y*y)* TMath::Exp(em*alxm - TMath::LnGamma(em + 1.0) - g); 00445 } while( Rndm() > t ); 00446 00447 return em; 00448 00449 } else { 00450 // use Gaussian approximation vor very large values 00451 return Gaus(0,1)*TMath::Sqrt(mean) + mean +0.5; 00452 } 00453 } 00454 00455 //______________________________________________________________________________ 00456 void TRandom::Rannor(Float_t &a, Float_t &b) 00457 { 00458 // Return 2 numbers distributed following a gaussian with mean=0 and sigma=1 00459 00460 Double_t r, x, y, z; 00461 00462 y = Rndm(); 00463 z = Rndm(); 00464 x = z * 6.28318530717958623; 00465 r = TMath::Sqrt(-2*TMath::Log(y)); 00466 a = (Float_t)(r * TMath::Sin(x)); 00467 b = (Float_t)(r * TMath::Cos(x)); 00468 } 00469 00470 //______________________________________________________________________________ 00471 void TRandom::Rannor(Double_t &a, Double_t &b) 00472 { 00473 // Return 2 numbers distributed following a gaussian with mean=0 and sigma=1 00474 00475 Double_t r, x, y, z; 00476 00477 y = Rndm(); 00478 z = Rndm(); 00479 x = z * 6.28318530717958623; 00480 r = TMath::Sqrt(-2*TMath::Log(y)); 00481 a = r * TMath::Sin(x); 00482 b = r * TMath::Cos(x); 00483 } 00484 00485 //_____________________________________________________________________________ 00486 void TRandom::ReadRandom(const char *filename) 00487 { 00488 // 00489 // Reads saved random generator status from filename 00490 // 00491 if (!gDirectory) return; 00492 char *fntmp = gSystem->ExpandPathName(filename); 00493 TDirectory *file = (TDirectory*)gROOT->ProcessLine(Form("TFile::Open(\"%s\");",fntmp)); 00494 delete [] fntmp; 00495 if(file && file->GetFile()) { 00496 gDirectory->ReadTObject(this,GetName()); 00497 delete file; 00498 } 00499 } 00500 00501 //______________________________________________________________________________ 00502 Double_t TRandom::Rndm(Int_t) 00503 { 00504 // Machine independent random number generator. 00505 // Based on the BSD Unix (Rand) Linear congrential generator 00506 // Produces uniformly-distributed floating points between 0 and 1. 00507 // Identical sequence on all machines of >= 32 bits. 00508 // Periodicity = 2**31 00509 // generates a number in ]0,1] 00510 // Note that this is a generator which is known to have defects 00511 // (the lower random bits are correlated) and therefore should NOT be 00512 // used in any statistical study. 00513 00514 #ifdef OLD_TRANDOM_IMPL 00515 const Double_t kCONS = 4.6566128730774E-10; 00516 const Int_t kMASK24 = 2147483392; 00517 00518 fSeed *= 69069; 00519 UInt_t jy = (fSeed&kMASK24); // Set lower 8 bits to zero to assure exact float 00520 if (jy) return kCONS*jy; 00521 return Rndm(); 00522 #endif 00523 00524 const Double_t kCONS = 4.6566128730774E-10; // (1/pow(2,31)) 00525 fSeed = (1103515245 * fSeed + 12345) & 0x7fffffffUL; 00526 00527 if (fSeed) return kCONS*fSeed; 00528 return Rndm(); 00529 } 00530 00531 //______________________________________________________________________________ 00532 void TRandom::RndmArray(Int_t n, Double_t *array) 00533 { 00534 // Return an array of n random numbers uniformly distributed in ]0,1] 00535 00536 const Double_t kCONS = 4.6566128730774E-10; // (1/pow(2,31)) 00537 Int_t i=0; 00538 while (i<n) { 00539 fSeed = (1103515245 * fSeed + 12345) & 0x7fffffffUL; 00540 if (fSeed) {array[i] = kCONS*fSeed; i++;} 00541 } 00542 } 00543 00544 //______________________________________________________________________________ 00545 void TRandom::RndmArray(Int_t n, Float_t *array) 00546 { 00547 // Return an array of n random numbers uniformly distributed in ]0,1] 00548 00549 const Double_t kCONS = 4.6566128730774E-10; // (1/pow(2,31)) 00550 const Int_t kMASK24 = 0x7fffff00; 00551 UInt_t jy; 00552 Int_t i=0; 00553 while (i<n) { 00554 fSeed = (1103515245 * fSeed + 12345) & 0x7fffffffUL; 00555 jy = (fSeed&kMASK24); // Set lower 8 bits to zero to assure exact float 00556 if (fSeed) {array[i] = Float_t(kCONS*fSeed); i++;} 00557 } 00558 } 00559 00560 //______________________________________________________________________________ 00561 void TRandom::SetSeed(UInt_t seed) 00562 { 00563 // Set the random generator seed. Note that default value is zero, which is different than the 00564 // default value used when constructing the class. 00565 // If the seed is zero the seed is set to a random value 00566 // which in case of TRandom depends on the machine clock. 00567 // Note that the machine clock is returned with a precision of 1 second. 00568 // If one calls SetSeed(0) within a loop and the loop time is less than 1s, 00569 // all generated numbers will be identical! 00570 // Instead if a different generator implementation is used (TRandom1 , 2 or 3) the seed is generated using 00571 // a 128 bit UUID. This results in different seeds and then random sequence for every SetSeed(0) call. 00572 00573 if( seed==0 ) { 00574 time_t curtime; // Set 'random' seed number if seed=0 00575 time(&curtime); // Get current time in fSeed. 00576 fSeed = (UInt_t)curtime; 00577 } else { 00578 fSeed = seed; 00579 } 00580 } 00581 00582 //______________________________________________________________________________ 00583 void TRandom::Sphere(Double_t &x, Double_t &y, Double_t &z, Double_t r) 00584 { 00585 // generates random vectors, uniformly distributed over the surface 00586 // of a sphere of given radius. 00587 // Input : r = sphere radius 00588 // Output: x,y,z a random 3-d vector of length r 00589 // Method: (based on algorithm suggested by Knuth and attributed to Robert E Knop) 00590 // which uses less random numbers than the CERNLIB RN23DIM algorithm 00591 00592 Double_t a=0,b=0,r2=1; 00593 while (r2 > 0.25) { 00594 a = Rndm() - 0.5; 00595 b = Rndm() - 0.5; 00596 r2 = a*a + b*b; 00597 } 00598 z = r* ( -1. + 8.0 * r2 ); 00599 00600 Double_t scale = 8.0 * r * TMath::Sqrt(0.25 - r2); 00601 x = a*scale; 00602 y = b*scale; 00603 } 00604 00605 //______________________________________________________________________________ 00606 Double_t TRandom::Uniform(Double_t x1) 00607 { 00608 // returns a uniform deviate on the interval ]0, x1]. 00609 00610 Double_t ans = Rndm(); 00611 return x1*ans; 00612 } 00613 00614 //______________________________________________________________________________ 00615 Double_t TRandom::Uniform(Double_t x1, Double_t x2) 00616 { 00617 // returns a uniform deviate on the interval ]x1, x2]. 00618 00619 Double_t ans= Rndm(); 00620 return x1 + (x2-x1)*ans; 00621 } 00622 00623 //_____________________________________________________________________________ 00624 void TRandom::WriteRandom(const char *filename) 00625 { 00626 // 00627 // Writes random generator status to filename 00628 // 00629 if (!gDirectory) return; 00630 char *fntmp = gSystem->ExpandPathName(filename); 00631 TDirectory *file = (TDirectory*)gROOT->ProcessLine(Form("TFile::Open(\"%s\",\"recreate\");",fntmp)); 00632 delete [] fntmp; 00633 if(file && file->GetFile()) { 00634 gDirectory->WriteTObject(this,GetName()); 00635 delete file; 00636 } 00637 }