ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TGenPhaseSpace.cxx
Go to the documentation of this file.
1 // @(#)root/physics:$Id$
2 // Author: Rene Brun , Valerio Filippini 06/09/2000
3 
4 //_____________________________________________________________________________________
5 //
6 // Utility class to generate n-body event,
7 // with constant cross-section (default)
8 // or with Fermi energy dependence (opt="Fermi").
9 // The event is generated in the center-of-mass frame,
10 // but the decay products are finally boosted
11 // using the betas of the original particle.
12 //
13 // The code is based on the GENBOD function (W515 from CERNLIB)
14 // using the Raubold and Lynch method
15 // F. James, Monte Carlo Phase Space, CERN 68-15 (1968)
16 //
17 // see example of use in $ROOTSYS/tutorials/physics/PhaseSpace.C
18 //
19 // Note that Momentum, Energy units are Gev/C, GeV
20 
21 #include "TGenPhaseSpace.h"
22 #include "TRandom.h"
23 #include "TMath.h"
24 #include <stdlib.h>
25 
26 const Int_t kMAXP = 18;
27 
29 
30 ////////////////////////////////////////////////////////////////////////////////
31 ///the PDK function
32 
34 {
35  Double_t x = (a-b-c)*(a+b+c)*(a-b+c)*(a+b-c);
36  x = TMath::Sqrt(x)/(2*a);
37  return x;
38 }
39 
40 ////////////////////////////////////////////////////////////////////////////////
41 ///special max function
42 
43 Int_t DoubleMax(const void *a, const void *b)
44 {
45  Double_t aa = * ((Double_t *) a);
46  Double_t bb = * ((Double_t *) b);
47  if (aa > bb) return 1;
48  if (aa < bb) return -1;
49  return 0;
50 
51 }
52 
53 ////////////////////////////////////////////////////////////////////////////////
54 ///copy constructor
55 
57 {
58  fNt = gen.fNt;
59  fWtMax = gen.fWtMax;
60  fTeCmTm = gen.fTeCmTm;
61  fBeta[0] = gen.fBeta[0];
62  fBeta[1] = gen.fBeta[1];
63  fBeta[2] = gen.fBeta[2];
64  for (Int_t i=0;i<fNt;i++) {
65  fMass[i] = gen.fMass[i];
66  fDecPro[i] = gen.fDecPro[i];
67  }
68 }
69 
70 
71 ////////////////////////////////////////////////////////////////////////////////
72 /// Assignment operator
73 
75 {
76  TObject::operator=(gen);
77  fNt = gen.fNt;
78  fWtMax = gen.fWtMax;
79  fTeCmTm = gen.fTeCmTm;
80  fBeta[0] = gen.fBeta[0];
81  fBeta[1] = gen.fBeta[1];
82  fBeta[2] = gen.fBeta[2];
83  for (Int_t i=0;i<fNt;i++) {
84  fMass[i] = gen.fMass[i];
85  fDecPro[i] = gen.fDecPro[i];
86  }
87  return *this;
88 }
89 
90 ////////////////////////////////////////////////////////////////////////////////
91 /// Generate a random final state.
92 /// The function returns the weigth of the current event.
93 /// The TLorentzVector of each decay product can be obtained using GetDecay(n).
94 ///
95 /// Note that Momentum, Energy units are Gev/C, GeV
96 
98 {
99  Double_t rno[kMAXP];
100  rno[0] = 0;
101  Int_t n;
102  if (fNt>2) {
103  for (n=1; n<fNt-1; n++) rno[n]=gRandom->Rndm(); // fNt-2 random numbers
104  qsort(rno+1 ,fNt-2 ,sizeof(Double_t) ,DoubleMax); // sort them
105  }
106  rno[fNt-1] = 1;
107 
108  Double_t invMas[kMAXP], sum=0;
109  for (n=0; n<fNt; n++) {
110  sum += fMass[n];
111  invMas[n] = rno[n]*fTeCmTm + sum;
112  }
113 
114  //
115  //-----> compute the weight of the current event
116  //
117  Double_t wt=fWtMax;
118  Double_t pd[kMAXP];
119  for (n=0; n<fNt-1; n++) {
120  pd[n] = PDK(invMas[n+1],invMas[n],fMass[n+1]);
121  wt *= pd[n];
122  }
123 
124  //
125  //-----> complete specification of event (Raubold-Lynch method)
126  //
127  fDecPro[0].SetPxPyPzE(0, pd[0], 0 , TMath::Sqrt(pd[0]*pd[0]+fMass[0]*fMass[0]) );
128 
129  Int_t i=1;
130  Int_t j;
131  while (1) {
132  fDecPro[i].SetPxPyPzE(0, -pd[i-1], 0 , TMath::Sqrt(pd[i-1]*pd[i-1]+fMass[i]*fMass[i]) );
133 
134  Double_t cZ = 2*gRandom->Rndm() - 1;
135  Double_t sZ = TMath::Sqrt(1-cZ*cZ);
136  Double_t angY = 2*TMath::Pi() * gRandom->Rndm();
137  Double_t cY = TMath::Cos(angY);
138  Double_t sY = TMath::Sin(angY);
139  for (j=0; j<=i; j++) {
140  TLorentzVector *v = fDecPro+j;
141  Double_t x = v->Px();
142  Double_t y = v->Py();
143  v->SetPx( cZ*x - sZ*y );
144  v->SetPy( sZ*x + cZ*y ); // rotation around Z
145  x = v->Px();
146  Double_t z = v->Pz();
147  v->SetPx( cY*x - sY*z );
148  v->SetPz( sY*x + cY*z ); // rotation around Y
149  }
150 
151  if (i == (fNt-1)) break;
152 
153  Double_t beta = pd[i] / sqrt(pd[i]*pd[i] + invMas[i]*invMas[i]);
154  for (j=0; j<=i; j++) fDecPro[j].Boost(0,beta,0);
155  i++;
156  }
157 
158  //
159  //---> final boost of all particles
160  //
161  for (n=0;n<fNt;n++) fDecPro[n].Boost(fBeta[0],fBeta[1],fBeta[2]);
162 
163  //
164  //---> return the weigth of event
165  //
166  return wt;
167 }
168 
169 ////////////////////////////////////////////////////////////////////////////////
170 ///return Lorentz vector corresponding to decay n
171 
173 {
174  if (n>fNt) return 0;
175  return fDecPro+n;
176 }
177 
178 
179 ////////////////////////////////////////////////////////////////////////////////
180 /// input:
181 /// TLorentzVector &P: decay particle (Momentum, Energy units are Gev/C, GeV)
182 /// Int_t nt: number of decay products
183 /// Double_t *mass: array of decay product masses
184 /// Option_t *opt: default -> constant cross section
185 /// "Fermi" -> Fermi energy dependece
186 /// return value:
187 /// kTRUE: the decay is permitted by kinematics
188 /// kFALSE: the decay is forbidden by kinematics
189 ///
190 
192  const Double_t *mass, Option_t *opt)
193 {
194  Int_t n;
195  fNt = nt;
196  if (fNt<2 || fNt>18) return kFALSE; // no more then 18 particle
197 
198  //
199  //
200  //
201  fTeCmTm = P.Mag(); // total energy in C.M. minus the sum of the masses
202  for (n=0;n<fNt;n++) {
203  fMass[n] = mass[n];
204  fTeCmTm -= mass[n];
205  }
206 
207  if (fTeCmTm<=0) return kFALSE; // not enough energy for this decay
208 
209  //
210  //------> the max weigth depends on opt:
211  // opt == "Fermi" --> fermi energy dependence for cross section
212  // else --> constant cross section as function of TECM (default)
213  //
214  if (strcasecmp(opt,"fermi")==0) {
215  // ffq[] = pi * (2*pi)**(FNt-2) / (FNt-2)!
216  Double_t ffq[] = {0
217  ,3.141592, 19.73921, 62.01255, 129.8788, 204.0131
218  ,256.3704, 268.4705, 240.9780, 189.2637
219  ,132.1308, 83.0202, 47.4210, 24.8295
220  ,12.0006, 5.3858, 2.2560, 0.8859 };
221  fWtMax = TMath::Power(fTeCmTm,fNt-2) * ffq[fNt-1] / P.Mag();
222 
223  } else {
224  Double_t emmax = fTeCmTm + fMass[0];
225  Double_t emmin = 0;
226  Double_t wtmax = 1;
227  for (n=1; n<fNt; n++) {
228  emmin += fMass[n-1];
229  emmax += fMass[n];
230  wtmax *= PDK(emmax, emmin, fMass[n]);
231  }
232  fWtMax = 1/wtmax;
233  }
234 
235  //
236  //----> save the betas of the decaying particle
237  //
238  if (P.Beta()) {
239  Double_t w = P.Beta()/P.Rho();
240  fBeta[0] = P(0)*w;
241  fBeta[1] = P(1)*w;
242  fBeta[2] = P(2)*w;
243  }
244  else fBeta[0]=fBeta[1]=fBeta[2]=0;
245 
246  return kTRUE;
247 }
void SetPxPyPzE(Double_t px, Double_t py, Double_t pz, Double_t e)
void SetPx(Double_t a)
Double_t fMass[18]
virtual Double_t Rndm(Int_t i=0)
Machine independent random number generator.
Definition: TRandom.cxx:512
return c
const char Option_t
Definition: RtypesCore.h:62
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
const Bool_t kFALSE
Definition: Rtypes.h:92
const Int_t kMAXP
Double_t Generate()
Generate a random final state.
ClassImp(TGenPhaseSpace) Double_t TGenPhaseSpace
the PDK function
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:501
double beta(double x, double y)
Calculates the beta function.
double sqrt(double)
Double_t x[n]
Definition: legend1.C:17
TGenPhaseSpace & operator=(const TGenPhaseSpace &gen)
Assignment operator.
TLorentzVector fDecPro[18]
TObject & operator=(const TObject &rhs)
TObject assignment operator.
Definition: TObject.cxx:102
void SetPy(Double_t a)
Int_t DoubleMax(const void *a, const void *b)
special max function
Float_t z[5]
Definition: Ifit.C:16
Double_t Rho() const
SVector< double, 2 > v
Definition: Dict.h:5
tuple w
Definition: qtexample.py:51
R__EXTERN TRandom * gRandom
Definition: TRandom.h:62
void SetPz(Double_t a)
Double_t Cos(Double_t)
Definition: TMath.h:424
Double_t Pi()
Definition: TMath.h:44
Double_t Px() const
Double_t Py() const
MyComplex< T > P(MyComplex< T > z, T c_real, T c_imag)
[MyComplex]
Definition: mandel.cpp:155
double Double_t
Definition: RtypesCore.h:55
Bool_t SetDecay(TLorentzVector &P, Int_t nt, const Double_t *mass, Option_t *opt="")
input: TLorentzVector &P: decay particle (Momentum, Energy units are Gev/C, GeV) Int_t nt: number of ...
Double_t y[n]
Definition: legend1.C:17
Mother of all ROOT objects.
Definition: TObject.h:58
Double_t Mag() const
Double_t PDK(Double_t a, Double_t b, Double_t c)
TLorentzVector * GetDecay(Int_t n)
return Lorentz vector corresponding to decay n
Double_t Beta() const
Double_t Sin(Double_t)
Definition: TMath.h:421
Double_t fTeCmTm
Double_t Pz() const
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
const Bool_t kTRUE
Definition: Rtypes.h:91
const Int_t n
Definition: legend1.C:16
Double_t fBeta[3]
TText gen(9, 20.7,"rootcint")