library: libCore #include "TRandom.h" |
TRandom
class description - header file - source file - inheritance tree (.pdf)
public:
TRandom(UInt_t seed = 65539)
TRandom(const TRandom&)
virtual ~TRandom()
virtual Int_t Binomial(Int_t ntot, Double_t prob)
virtual Double_t BreitWigner(Double_t mean = 0, Double_t gamma = 1)
virtual void Circle(Double_t& x, Double_t& y, Double_t r)
static TClass* Class()
virtual Double_t Exp(Double_t tau)
virtual Double_t Gaus(Double_t mean = 0, Double_t sigma = 1)
virtual UInt_t GetSeed()
virtual UInt_t Integer(UInt_t imax)
virtual TClass* IsA() const
virtual Double_t Landau(Double_t mean = 0, Double_t sigma = 1)
TRandom& operator=(const TRandom&)
virtual Int_t Poisson(Double_t mean)
virtual Double_t PoissonD(Double_t mean)
virtual void Rannor(Float_t& a, Float_t& b)
virtual void Rannor(Double_t& a, Double_t& b)
virtual void ReadRandom(const char* filename)
virtual Double_t Rndm(Int_t i = 0)
virtual void RndmArray(Int_t n, Float_t* array)
virtual void RndmArray(Int_t n, Double_t* array)
virtual void SetSeed(UInt_t seed = 65539)
virtual void ShowMembers(TMemberInspector& insp, char* parent)
virtual void Sphere(Double_t& x, Double_t& y, Double_t& z, Double_t r)
virtual void Streamer(TBuffer& b)
void StreamerNVirtual(TBuffer& b)
virtual Double_t Uniform(Double_t x1 = 1)
virtual Double_t Uniform(Double_t x1, Double_t x2)
virtual void WriteRandom(const char* filename)
protected:
UInt_t fSeed Random number generator seed
TRandom
basic Random number generator class (periodicity = 10**9).
Note that this is a very simple generator (linear congruential)
which is known to have defects (the lower random bits are correlated)
and therefore should NOT be used in any statistical study.
One should use instead TRandom1, TRandom2 or TRandom3.
TRandom3, is based on the "Mersenne Twister generator", and is the recommended one,
since it has good random proprieties (period of about 10**6000 ) and it is fast.
TRandom1, based on the RANLUX algorithm, has mathematically proven random proprieties
and a period of about 10**171. It is however slower than the others.
TRandom2, is based on the Tausworthe generator of L'Ecuyer, and it has the advantage
of being fast and using only 3 words (of 32 bits) for the state. The period is 10**26.
The following table shows some timings (in nanoseconds/call)
for the random numbers obtained using an Intel Pentium 3.0 GHz running Linux
and using the gcc 3.2.3 compiler
TRandom 34 ns/call (BAD Generator)
TRandom1 242 ns/call
TRandom2 37 ns/call
TRandom3 45 ns/call
The following basic Random distributions are provided:
===================================================
-Exp(tau)
-Integer(imax)
-Gaus(mean,sigma)
-Rndm()
-Uniform(x1)
-Landau(mpv,sigma)
-Poisson(mean)
-Binomial(ntot,prob)
Random numbers distributed according to 1-d, 2-d or 3-d distributions
=====================================================================
contained in TF1, TF2 or TF3 objects.
For example, to get a random number distributed following abs(sin(x)/x)*sqrt(x)
you can do:
TF1 *f1 = new TF1("f1","abs(sin(x)/x)*sqrt(x)",0,10);
double r = f1->GetRandom();
The technique of using a TF1,2 or 3 function is very powerful.
It is also more precise than using the basic functions (except Rndm).
With a TF1 function, for example, the real integral of the function
is correctly calculated in the specified range of the function.
Getting a number from a TF1 function is also very fast.
The following table shows some timings (in microsecons/call)
for basic functions and TF1 functions.
The left column is with the compiler, the right column with CINT.
Numbers have been obtained on a Pentium 233Mhz running Linux.
g++ CINT
Rndm.............. 0.330 4.15
Gaus.............. 2.220 6.77
Landau............ 21.590 46.82
Binomial(5,0.5)... 0.890 5.34
Binomial(15,0.5).. 0.920 5.36
Poisson(3)........ 2.170 5.93
Poisson(10)....... 4.160 7.95
Poisson(70)....... 21.510 25.27
Poisson(100)...... 2.910 6.72
GausTF1........... 2.070 4.73
LandauTF1......... 2.100 4.73
Note that the time to generate a number from an arbitrary TF1 function
is independent of the complexity of the function.
For Landau distribution, it is recommended to use the TF1 technique.
TH1::FillRandom(TH1 *) or TH1::FillRandom(const char *tf1name)
==============================================================
can be used to fill an histogram (1-d, 2-d, 3-d from an existing histogram
or from an existing function.
Note this interesting feature when working with objects
=======================================================
You can use several TRandom objects, each with their "independent"
random sequence. For example, one can imagine
TRandom *eventGenerator = new TRandom();
TRandom *tracking = new TRandom();
eventGenerator can be used to generate the event kinematics.
tracking can be used to track the generated particles with random numbers
independent from eventGenerator.
This very interesting feature gives the possibility to work with simple
and very fast random number generators without worrying about
random number periodicity as it was the case with Fortran.
One can use TRandom::SetSeed to modify the seed of one generator.
a TRandom object may be written to a Root file
==============================================
-as part of another object
-or with its own key (example gRandom->Write("Random");
The small program below has been used to get the values in the table above.
#ifndef __CINT__
#include "TROOT.h"
#include "TF1.h"
#include "TRandom.h"
#include "TStopwatch.h"
void rand();
//______________________________________________________________________________
int main()
{
TROOT simple("simple","Test of random numbers");
rand();
}
#endif
void rand() {
int i, N = 1000000;
double cpn = 1000000./N;
double x;
TStopwatch sw;
sw.Start();
for (i=0;i<N;i++) {
x = gRandom->Rndm(i);
}
printf("Rndm.............. %8.3f microseconds/call\n",sw.CpuTime()*cpn);
sw.Start();
for (i=0;i<N;i++) {
x = gRandom->Gaus(0,1);
}
printf("Gaus.............. %8.3f\n",sw.CpuTime()*cpn);
sw.Start();
for (i=0;i<N;i++) {
x = gRandom->Landau(0,1);
}
printf("Landau............ %8.3f\n",sw.CpuTime()*cpn);
sw.Start();
for (i=0;i<N;i++) {
x = gRandom->Binomial(5,0.5);
}
printf("Binomial(5,0.5)... %8.3f\n",sw.CpuTime()*cpn);
sw.Start();
for (i=0;i<N;i++) {
x = gRandom->Binomial(15,0.5);
}
printf("Binomial(15,0.5).. %8.3f\n",sw.CpuTime()*cpn);
sw.Start();
for (i=0;i<N;i++) {
x = gRandom->Poisson(3);
}
printf("Poisson(3)........ %8.3f\n",sw.CpuTime()*cpn);
sw.Start();
for (i=0;i<N;i++) {
x = gRandom->Poisson(10);
}
printf("Poisson(10)....... %8.3f\n",sw.CpuTime()*cpn);
sw.Start();
for (i=0;i<N;i++) {
x = gRandom->Poisson(70);
}
printf("Poisson(70)....... %8.3f\n",sw.CpuTime()*cpn);
sw.Start();
for (i=0;i<N;i++) {
x = gRandom->Poisson(100);
}
printf("Poisson(100)...... %8.3f\n",sw.CpuTime()*cpn);
TF1 *f1 = new TF1("f1","gaus",-4,4);
f1->SetParameters(1,0,1);
sw.Start();
for (i=0;i<N;i++) {
x = f1->GetRandom();
}
printf("GausTF1........... %8.3f\n",sw.CpuTime()*cpn);
TF1 *f2 = new TF1("f2","landau",-5,15);
f2->SetParameters(1,0,1);
sw.Start();
for (i=0;i<N;i++) {
x = f2->GetRandom();
}
printf("LandauTF1......... %8.3f\n",sw.CpuTime()*cpn);
}
TRandom(UInt_t seed)
*-*-*-*-*-*-*-*-*-*-*default constructor*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
*-* ===================
~TRandom()
*-*-*-*-*-*-*-*-*-*-*default destructor*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
*-* ==================
Int_t Binomial(Int_t ntot, Double_t prob)
Generates a random integer N according to the binomial law
Coded from Los Alamos report LA-5061-MS
N is binomially distributed between 0 and ntot inclusive
with mean prob*ntot.
prob is between 0 and 1.
Note: This function should not be used when ntot is large (say >100).
The normal approximation is then recommended instead
(with mean =*ntot+0.5 and standard deviation sqrt(ntot*prob*(1-prob)).
void Circle(Double_t &x, Double_t &y, Double_t r)
generates random vectors, uniformly distributed over a circle of given radius.
Input : r = circle radius
Output: x,y a random 2-d vector of length r
Double_t Landau(Double_t mpv, Double_t sigma)
Generate a random number following a Landau distribution
with mpv(most probable value) and sigma
Converted by Rene Brun from CERNLIB routine ranlan(G110)
Int_t Poisson(Double_t mean)
Generates a random integer N according to a Poisson law.
Prob(N) = exp(-mean)*mean^N/Factorial(N)
Use a different procedure according to the mean value.
The algorithm is the same used by CLHEP
For lower value (mean < 25) use the rejection method based on
the exponential
For higher values use a rejection method comparing with a Lorentzian
distribution, as suggested by several authors
This routine since is returning 32 bits integer will not work for values larger than 2*10**9
One should then use the Trandom::PoissonD for such large values
void Rannor(Float_t &a, Float_t &b)
Return 2 numbers distributed following a gaussian with mean=0 and sigma=1
void Rannor(Double_t &a, Double_t &b)
Return 2 numbers distributed following a gaussian with mean=0 and sigma=1
void ReadRandom(const char *filename)
Reads saved random generator status from filename
Double_t Rndm(Int_t)
Machine independent random number generator.
Based on the BSD Unix (Rand) Linear congrential generator
Produces uniformly-distributed floating points between 0 and 1.
Identical sequence on all machines of >= 32 bits.
Periodicity = 2**31
generates a number in ]0,1]
Note that this is a generator which is known to have defects
(the lower random bits are correlated) and therefore should NOT be
used in any statistical study.
void RndmArray(Int_t n, Float_t *array)
Return an array of n random numbers uniformly distributed in ]0,1]
void SetSeed(UInt_t seed)
Set the random generator seed
if seed is zero, the seed is set to the current machine clock
Note that the machine clock is returned with a precision of 1 second.
If one calls SetSeed(0) within a loop and the loop time is less than 1s,
all generated numbers will be identical!
void Sphere(Double_t &x, Double_t &y, Double_t &z, Double_t r)
generates random vectors, uniformly distributed over the surface
of a sphere of given radius.
Input : r = sphere radius
Output: x,y,z a random 3-d vector of length r
Method: (based on algorithm suggested by Knuth and attributed to Robert E Knop)
which uses less random numbers than the CERNLIB RN23DIM algorithm
void WriteRandom(const char *filename)
Writes random generator status to filename
Author: Rene Brun, Lorenzo Moneta 15/12/95
Last update: root/base:$Name: $:$Id: TRandom.cxx,v 1.33 2006/05/28 06:15:54 brun Exp $
Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
ROOT page - Class index - Class Hierarchy - Top of the page
This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.