TGenPhaseSpace.cxx

Go to the documentation of this file.
00001 // @(#)root/physics:$Id: TGenPhaseSpace.cxx 20882 2007-11-19 11:31:26Z rdm $
00002 // Author: Rene Brun , Valerio Filippini  06/09/2000 
00003 
00004 //_____________________________________________________________________________________
00005 //
00006 //  Utility class to generate n-body event,
00007 //  with constant cross-section (default)
00008 //  or with Fermi energy dependence (opt="Fermi").
00009 //  The event is generated in the center-of-mass frame, 
00010 //  but the decay products are finally boosted
00011 //  using the betas of the original particle.
00012 //
00013 //  The code is based on the GENBOD function (W515 from CERNLIB)
00014 //  using the Raubold and Lynch method
00015 //      F. James, Monte Carlo Phase Space, CERN 68-15 (1968)
00016 //
00017 // see example of use in $ROOTSYS/tutorials/physics/PhaseSpace.C
00018 //
00019 // Note that Momentum, Energy units are Gev/C, GeV
00020 
00021 #include "TGenPhaseSpace.h"
00022 #include "TRandom.h"
00023 #include "TMath.h"
00024 #include <stdlib.h>
00025 
00026 const Int_t kMAXP = 18;
00027 
00028 ClassImp(TGenPhaseSpace)
00029 
00030 //_____________________________________________________________________________________
00031 Double_t TGenPhaseSpace::PDK(Double_t a, Double_t b, Double_t c) 
00032 {
00033    //the PDK function
00034    Double_t x = (a-b-c)*(a+b+c)*(a-b+c)*(a+b-c);
00035    x = TMath::Sqrt(x)/(2*a);
00036    return x;
00037 }
00038 
00039 //_____________________________________________________________________________________
00040 Int_t DoubleMax(const void *a, const void *b) 
00041 {
00042    //special max function
00043    Double_t aa = * ((Double_t *) a); 
00044    Double_t bb = * ((Double_t *) b); 
00045    if (aa > bb) return  1;
00046    if (aa < bb) return -1;
00047    return 0;
00048 
00049 }
00050 
00051 //__________________________________________________________________________________________________
00052 TGenPhaseSpace::TGenPhaseSpace(const TGenPhaseSpace &gen) : TObject(gen)
00053 {
00054    //copy constructor
00055    fNt      = gen.fNt;
00056    fWtMax   = gen.fWtMax;
00057    fTeCmTm  = gen.fTeCmTm;
00058    fBeta[0] = gen.fBeta[0];
00059    fBeta[1] = gen.fBeta[1];
00060    fBeta[2] = gen.fBeta[2];
00061    for (Int_t i=0;i<fNt;i++) {
00062       fMass[i]   = gen.fMass[i];
00063       fDecPro[i] = gen.fDecPro[i];
00064    }
00065 }
00066 
00067 //__________________________________________________________________________________________________
00068 Double_t TGenPhaseSpace::Generate() 
00069 {
00070    //  Generate a random final state.
00071    //  The function returns the weigth of the current event.
00072    //  The TLorentzVector of each decay product can be obtained using GetDecay(n).
00073    //
00074    // Note that Momentum, Energy units are Gev/C, GeV
00075 
00076    Double_t rno[kMAXP];
00077    rno[0] = 0;
00078    Int_t n;
00079    if (fNt>2) {
00080       for (n=1; n<fNt-1; n++)  rno[n]=gRandom->Rndm();   // fNt-2 random numbers
00081       qsort(rno+1 ,fNt-2 ,sizeof(Double_t) ,DoubleMax);  // sort them
00082    }
00083    rno[fNt-1] = 1;
00084 
00085    Double_t invMas[kMAXP], sum=0;
00086    for (n=0; n<fNt; n++) {
00087       sum      += fMass[n];
00088       invMas[n] = rno[n]*fTeCmTm + sum;
00089    }
00090 
00091    //
00092    //-----> compute the weight of the current event
00093    //
00094    Double_t wt=fWtMax;
00095    Double_t pd[kMAXP];
00096    for (n=0; n<fNt-1; n++) {
00097       pd[n] = PDK(invMas[n+1],invMas[n],fMass[n+1]);
00098       wt *= pd[n];
00099    }
00100 
00101    //
00102    //-----> complete specification of event (Raubold-Lynch method)
00103    //
00104    fDecPro[0].SetPxPyPzE(0, pd[0], 0 , TMath::Sqrt(pd[0]*pd[0]+fMass[0]*fMass[0]) );
00105 
00106    Int_t i=1;
00107    Int_t j;
00108    while (1) {
00109       fDecPro[i].SetPxPyPzE(0, -pd[i-1], 0 , TMath::Sqrt(pd[i-1]*pd[i-1]+fMass[i]*fMass[i]) );
00110 
00111       Double_t cZ   = 2*gRandom->Rndm() - 1;
00112       Double_t sZ   = TMath::Sqrt(1-cZ*cZ);
00113       Double_t angY = 2*TMath::Pi() * gRandom->Rndm();
00114       Double_t cY   = TMath::Cos(angY);
00115       Double_t sY   = TMath::Sin(angY);
00116       for (j=0; j<=i; j++) {
00117          TLorentzVector *v = fDecPro+j;
00118          Double_t x = v->Px();
00119          Double_t y = v->Py();
00120          v->SetPx( cZ*x - sZ*y );
00121          v->SetPy( sZ*x + cZ*y );   // rotation around Z
00122          x = v->Px();
00123          Double_t z = v->Pz();
00124          v->SetPx( cY*x - sY*z );
00125          v->SetPz( sY*x + cY*z );   // rotation around Y
00126       }
00127 
00128       if (i == (fNt-1)) break;
00129 
00130       Double_t beta = pd[i] / sqrt(pd[i]*pd[i] + invMas[i]*invMas[i]);
00131       for (j=0; j<=i; j++) fDecPro[j].Boost(0,beta,0);
00132       i++;
00133    }
00134  
00135    //
00136    //---> final boost of all particles  
00137    //
00138    for (n=0;n<fNt;n++) fDecPro[n].Boost(fBeta[0],fBeta[1],fBeta[2]);
00139 
00140    //
00141    //---> return the weigth of event
00142    //
00143    return wt;
00144 }
00145 
00146 //__________________________________________________________________________________
00147 TLorentzVector *TGenPhaseSpace::GetDecay(Int_t n) 
00148 { 
00149    //return Lorentz vector corresponding to decay n
00150    if (n>fNt) return 0;
00151    return fDecPro+n;
00152 }
00153 
00154 
00155 //_____________________________________________________________________________________
00156 Bool_t TGenPhaseSpace::SetDecay(TLorentzVector &P, Int_t nt, 
00157    Double_t *mass, Option_t *opt) 
00158 {
00159    // input:
00160    // TLorentzVector &P:    decay particle (Momentum, Energy units are Gev/C, GeV)
00161    // Int_t nt:             number of decay products
00162    // Double_t *mass:       array of decay product masses
00163    // Option_t *opt:        default -> constant cross section
00164    //                       "Fermi" -> Fermi energy dependece
00165    // return value:
00166    // kTRUE:      the decay is permitted by kinematics
00167    // kFALSE:     the decay is forbidden by kinematics
00168    //
00169 
00170    Int_t n;
00171    fNt = nt;
00172    if (fNt<2 || fNt>18) return kFALSE;  // no more then 18 particle
00173 
00174    //
00175    //
00176    //
00177    fTeCmTm = P.Mag();           // total energy in C.M. minus the sum of the masses
00178    for (n=0;n<fNt;n++) {
00179       fMass[n]  = mass[n];
00180       fTeCmTm  -= mass[n];
00181    }
00182 
00183    if (fTeCmTm<=0) return kFALSE;    // not enough energy for this decay
00184 
00185    //
00186    //------> the max weigth depends on opt:
00187    //   opt == "Fermi"  --> fermi energy dependence for cross section
00188    //   else            --> constant cross section as function of TECM (default)
00189    //
00190    if (strcasecmp(opt,"fermi")==0) {  
00191       // ffq[] = pi * (2*pi)**(FNt-2) / (FNt-2)!
00192       Double_t ffq[] = {0 
00193                      ,3.141592, 19.73921, 62.01255, 129.8788, 204.0131
00194                      ,256.3704, 268.4705, 240.9780, 189.2637
00195                      ,132.1308,  83.0202,  47.4210,  24.8295
00196                      ,12.0006,   5.3858,   2.2560,   0.8859 };
00197       fWtMax = TMath::Power(fTeCmTm,fNt-2) * ffq[fNt-1] / P.Mag();
00198 
00199    } else {
00200       Double_t emmax = fTeCmTm + fMass[0];
00201       Double_t emmin = 0;
00202       Double_t wtmax = 1;
00203       for (n=1; n<fNt; n++) {
00204          emmin += fMass[n-1];
00205          emmax += fMass[n];
00206          wtmax *= PDK(emmax, emmin, fMass[n]);
00207       }
00208       fWtMax = 1/wtmax;
00209    }
00210 
00211    //
00212    //---->  save the betas of the decaying particle
00213    //
00214    if (P.Beta()) {
00215       Double_t w = P.Beta()/P.Rho();
00216       fBeta[0] = P(0)*w;
00217       fBeta[1] = P(1)*w;
00218       fBeta[2] = P(2)*w;
00219    }
00220    else fBeta[0]=fBeta[1]=fBeta[2]=0; 
00221 
00222    return kTRUE; 
00223 }

Generated on Tue Jul 5 14:37:46 2011 for ROOT_528-00b_version by  doxygen 1.5.1