]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenMUONLMR.cxx
Transport bit set.
[u/mrichter/AliRoot.git] / EVGEN / AliGenMUONLMR.cxx
1 //#include <TClonesArray.h> 
2
3 #include <TDatabasePDG.h>
4 #include <TFile.h>
5 #include "AliConst.h"
6 #include "AliGenMUONLMR.h" 
7 #include "AliMC.h" 
8 #include "AliRun.h" 
9 #include "AliGenEventHeader.h"
10
11 ClassImp(AliGenMUONLMR)
12
13
14 AliGenMUONLMR::AliGenMUONLMR () : AliGenMC(), fNMuMin(2), fGenSingleProc(-1), fCosTheta(0), fRhoLineShape(0), fHMultMu(0), fHNProc(0) { 
15   //
16   // default constructor 
17   //
18   // initialize pt and y distributions according to a fit to 
19   // Pythia simulation at sqrt(s) = 7 TeV
20   printf ("Using AliGenMUONLMR as generator\n");
21   for (Int_t ipart=0; ipart < fgkNpart; ipart++) fScaleMult[ipart] = 1; 
22   fScaleMult[kPionLMR] = 0; // set pion multiplicity to zero 
23   fScaleMult[kKaonLMR] = 0; // set kaon multiplicity to zero
24   Int_t pdg[7] = {211, 321, 221, 113, 223, 333, 331}; 
25   char* fptname[7] = {"fPtPion","fPtKaon","fPtEta","fPtRho","fPtOmega","fPtPhi","fPtEtaPrime"};
26   char* fyname[7] = {"fYPion","fYKaon","fYEta","fYRho","fYOmega","fYPhi","fYEtaPrime"}; 
27   char* fnname[7] = {"fMultPion","fMultKaon","fMultEta","fMultRho","fMultOmega","fMultPhi","fMultEtaPrime"};
28   char* fdname[2] = {"fDecPion","fDecKaon"};
29   Double_t ptparam[7][3] = { {1,0.427,2.52}, // pions from Pythia
30                              {1,0.58,2.57},  // kaons from Pythia
31                              {1,0.641,2.62},   // eta from Pythia
32                              {1,1.2,2.5},    // rho+omega from ALICE muon  
33                              {1,1.2,2.5},    // rho+omega from ALICE muon  
34                              {1,1.03,2.5},   // phi from ALICE muon  
35                              {1,0.72,2.5}};  // etaPrime from Pythia
36   
37   Double_t yparam[7][3] = { {1, 0.8251, 3.657}, // pions from pythia
38                             {1, 1.83, 2.698},   // kaons from pythia
39                             {1, 1.169, 3.282},  // eta from pythia
40                             {1, 1.234, 3.264},  // rho from pythia
41                             {1, 1.311, 3.223},  // omega from pythia
42                             {1, 2.388, 2.129},  // phi from pythia
43                             {1, 1.13,3.3}};     // eta prime from pythia
44
45   // multiplicity parameters from pythia
46   Double_t nparam[7][9] = { {353.582, 6.76263, 1.66979, 998.445, 9.73281, 12.6704, 175.187, 29.08, 40.2531},
47                             {1.e4, 0.2841,  0,0,0,0,0,0,0},
48                             {1.e4, 0.2647,  0,0,0,0,0,0,0},
49                             {7055, 0.1786,  0,0,0,0,0,0,0},
50                             {7500, 0.1896,  0,0,0,0,0,0,0},
51                             {5.e04, 1.167,  0,0,0,0,0,0,0}, 
52                             {2.9e04,0.714,  0,0,0,0,0,0,0}}; 
53
54   Double_t ctau[2] = {7.8045, 3.712};
55
56   for (Int_t i=0; i<fgkNpart; i++) { 
57     fPDG[i] = pdg[i]; 
58     if (i!=0) { 
59       fMult[i] = new TF1(fnname[i],"[0]*exp(-[1]*x)",0,30);
60       fMult[i]->SetParameters(nparam[i][0],nparam[i][1]);  
61     }
62     else { 
63       fMult[i] = new TF1(fnname[i],"gaus(0)+gaus(3)+gaus(6)",0,150);
64       for (Int_t j=0; j<9; j++) fMult[i]->SetParameter(j,nparam[i][j]);
65     }
66   
67     fPt[i] = new TF1(fptname[i],AliGenMUONLMR::PtDistr,0,20,3);
68     fPt[i]->SetParameters(ptparam[i][0], ptparam[i][1], ptparam[i][2]);  
69     fY[i] = new TF1(fyname[i],AliGenMUONLMR::YDistr,-10,10,3);
70     fY[i]->SetParameters(yparam[i][0], yparam[i][1], yparam[i][2]); 
71   }
72
73   for(Int_t i = 0; i<2; i++){
74     fDecay[i] = new TF1(fdname[i],"exp(-x/[0])",0,150);
75     fDecay[i]->SetParameter(0,ctau[i]);
76   }
77   
78   for (Int_t ipart = 0; ipart < fgkNpart; ipart++) { 
79     fParticle[ipart] = new TParticle(); 
80     fParticle[ipart]->SetPdgCode(fPDG[ipart]); 
81   }
82
83   TDatabasePDG *pdgdb = TDatabasePDG::Instance(); 
84   Double_t mumass = pdgdb->GetParticle(13)->Mass();
85   fMu[0] = new TParticle(); 
86   fMu[0]->SetPdgCode(-13); 
87   fMu[0]->SetCalcMass(mumass); 
88   fMu[1] = new TParticle(); 
89   fMu[1]->SetPdgCode(13); 
90   fMu[1]->SetCalcMass(mumass); 
91   
92   // function for polarized theta distributions
93   fCosTheta = new TF1 ("fCosTheta","1+[0]*x*x",-1,1);
94   fCosTheta->SetParameter(0,1);
95
96   // Dalitz decays 
97   Int_t nbins = 1000;
98   Double_t xmin = 0, xmax = 2; 
99   fDalitz[0] = new TH1F("hDalitzEta","",nbins,xmin,xmax);
100   fDalitz[1] = new TH1F("hDalitzOmega","",nbins,xmin,xmax);
101   fDalitz[2] = new TH1F("hDalitzEtaPrime","",nbins,xmin,xmax);
102   
103   Double_t meta   = pdgdb->GetParticle("eta")->Mass(); 
104   Double_t momega = pdgdb->GetParticle("omega")->Mass(); 
105   Double_t metaPrime = pdgdb->GetParticle("eta'")->Mass(); 
106   Double_t mpi0   = pdgdb->GetParticle("pi0")->Mass(); 
107   Double_t md3 = 0, mres = 0; 
108   
109   for (Int_t index = 0; index < 3; index++) { 
110     if (index == 0) { 
111       mres = meta; 
112       md3 = 0; 
113     }
114     else if (index == 1) { 
115       mres = momega; 
116       md3 = mpi0; 
117     }
118     else if (index == 2) { 
119       mres = metaPrime; 
120       md3 = 0; 
121     }
122     Double_t delta   = md3 * md3 / (mres * mres);
123     Double_t epsilon = mumass * mumass / (mres * mres);
124     Int_t nbins0 = fDalitz[index]->GetNbinsX();
125     Double_t xmin0 = fDalitz[index]->GetXaxis()->GetXmin(); 
126     Double_t deltax =  fDalitz[index]->GetBinWidth(1);
127     Double_t xd = xmin0 - deltax/2.; 
128     for (Int_t ibin = 0; ibin< nbins0; ibin++) { 
129       Double_t dalval = 0; 
130       xd += deltax; 
131       if (xd > 4. *epsilon) { 
132         Double_t bracket = TMath::Power(1. + xd/(1. - delta),2)      
133           - 4. * xd / ((1. - delta) * (1. - delta));
134         if (bracket > 0) { 
135           dalval = TMath::Power(bracket,1.5) /xd *
136             TMath::Sqrt(1 - 4 * epsilon / xd) * (1 + 2 * epsilon / xd) * 
137             FormFactor(xd * mres * mres, index);
138           fDalitz[index]->Fill(xd,dalval); 
139         }
140       }
141     }
142   }
143
144   fRhoLineShape = new TF1("fRhoLineShape",RhoLineShapeNew,0,2,2); 
145   fHMultMu = new TH1D("fHMultMu","Muon multiplicity",20,-0.5,19.5); 
146   fHNProc = new TH1D("fHNProc","Number of gen. evts. per process in 4 pi",9,-0.5,8.5); 
147 }
148
149 //-----------------------------------------------------------
150
151 AliGenMUONLMR::~AliGenMUONLMR()
152 {
153   // Default destructor
154   for (Int_t i=0; i<7; i++) { 
155     delete fPt[i]; 
156     delete fY[i]; 
157     delete fMult[i]; 
158     delete fParticle[i]; 
159   }    
160   
161   for (Int_t i=0; i<2; i++) { 
162     delete fDecay[i]; 
163     delete fMu[i]; 
164   }
165
166   for (Int_t i=0; i<3; i++) delete fDalitz[i]; 
167
168   delete fCosTheta; fCosTheta = 0;  
169   delete fRhoLineShape; fRhoLineShape = 0;  
170   delete fHMultMu; fHMultMu = 0;   
171   delete fHNProc;  fHNProc = 0;   
172 }
173
174 //-----------------------------------------------------------
175
176 void AliGenMUONLMR::FinishRun(){ 
177   // save some histograms to an output file 
178   Int_t nbins = fHNProc->GetNbinsX(); 
179   for (Int_t ibin=1; ibin <= nbins; ibin++) printf ("ibin = %d nEvProc = %g\n",
180                                                       ibin,fHNProc->GetBinContent(ibin));
181     TFile *fout = new TFile("AliGenMUONLMR_histos.root","recreate"); 
182   fHMultMu->Write(); 
183   fHNProc->Write(); 
184   fout->Close(); 
185 }
186
187 //-----------------------------------------------------------
188
189 Double_t AliGenMUONLMR::YDistr(Double_t *px, Double_t *par){ 
190   // function for rapidity distribution: plateau at par[0] +
191   // gaussian tails centered at par[1] and with par[2]=sigma  
192   Double_t x = TMath::Abs(px[0]);
193   Double_t func = 0;
194   if (x<par[1]) func = par[0]; 
195   else { 
196     Double_t z = (x-par[1])/(par[2]); 
197     func = par[0] * TMath::Exp(-0.5 * z * z); 
198   }
199   return func; 
200 }
201
202 //-----------------------------------------------------------
203
204 Double_t AliGenMUONLMR::PtDistr(Double_t *px, Double_t *par){
205   // pt distribution: power law 
206   Double_t x = px[0];
207   Double_t func = par[0] * x / TMath::Power((1+(x/par[1])*(x/par[1])),par[2]); 
208   return func; 
209 }
210
211 //-----------------------------------------------------------
212
213 void AliGenMUONLMR::Generate() {
214   //
215   // generate the low mass resonances and their decays according to  
216   // the multiplicity parameterized by pythia and BR from PDG  
217   // rapidity distributions parametrized from pythia 
218   // pt distributions from data (or pythia for etaprime) 
219   //
220   Double_t pxPushed[100], pyPushed[100], pzPushed[100], ePushed[100]; 
221   Int_t nmuons = -1, npartPushed = 0, pdgPushed[100]; 
222   Double_t polar[3]= {0,0,0};  // Polarisation of the parent particle (for GEANT tracking)
223   Double_t origin0[3];         // Origin of the generated parent particle (for GEANT tracking)
224   // Calculating vertex position per event
225   for (Int_t j=0;j<3;j++) origin0[j]=fOrigin[j];
226   if(fVertexSmear==kPerEvent) {
227     Vertex();
228     for (Int_t j=0;j<3;j++) origin0[j]=fVertex[j];
229   }
230   
231   printf ("In Generate()\n");
232   TParticle *mother; 
233   TDatabasePDG* pdg = TDatabasePDG::Instance();
234
235   Double_t pt, y, phi, mass, px, py, pz, ene, mt; 
236
237   const Int_t nproc = 9; 
238   Int_t idRes[nproc] = {kEtaLMR, kEtaLMR, kRhoLMR, kOmegaLMR, kOmegaLMR, kPhiLMR, kEtaPrimeLMR, kPionLMR, kKaonLMR}; 
239   Double_t BR[nproc] = {5.8e-6, 3.1e-4, 4.55e-5, 7.28e-4, 1.3e-4, 2.86e-4, 1.04e-4, 1, 0.6344};
240   //  Double_t BR[nproc] = {1, 1, 1, 1, 1, 1, 1, 1, 1};
241   Int_t idDec[nproc] = {0, 1, 0, 0, 1, 0, 1, 2, 2}; // 0:2body, 1:Dalitz, 2:pi/K 
242   Int_t mult[nproc] = {0,0,0,0,0,0,0,0,0}; 
243
244   while (nmuons < fNMuMin) { 
245
246     nmuons = 0; 
247     npartPushed = 0; 
248     for (Int_t iproc=0; iproc<nproc; iproc++) { 
249       if (fGenSingleProc == -1) { 
250         mult[iproc] = Int_t(fMult[idRes[iproc]]->GetRandom()*fScaleMult[idRes[iproc]]); 
251       }
252       else { 
253         if (iproc==fGenSingleProc) { 
254           mult[iproc] = 1; 
255           BR[iproc] = 1;
256         } 
257         else { 
258           mult[iproc] = 0; 
259           BR[iproc] = 0;
260         }
261       }
262     }
263     
264     if (fGenSingleProc == -1) { 
265       mult[1] = mult[0]; 
266       mult[4] = mult[3]; 
267     }
268     
269     for (Int_t iproc = 0; iproc < nproc; iproc++) { 
270 //       printf ("Multiplicity for process %d is %d\n",iproc,mult[iproc]); 
271       for (Int_t imult=0; imult<mult[iproc]; imult++) { 
272         if (gRandom->Rndm() < BR[iproc]) { 
273           fHNProc->Fill(iproc); 
274           Int_t ipart = idRes[iproc]; 
275           pt  = fPt[ipart]->GetRandom(); 
276           y   = fY[ipart]->GetRandom(); 
277           phi = gRandom->Rndm() * 2 * TMath::Pi(); 
278           mass = pdg->GetParticle(fPDG[ipart])->Mass(); 
279           px  = pt * TMath::Cos(phi); 
280           py  = pt * TMath::Sin(phi); 
281           mt  = TMath::Sqrt(pt * pt + mass * mass);
282           pz  = mt * TMath::SinH(y); 
283           ene = mt * TMath::CosH(y); 
284         
285           mother = fParticle[ipart]; 
286           mother->SetMomentum(px,py,pz,ene); 
287           mother->SetCalcMass(mass);
288           if (!KinematicSelection(mother,0)) continue; 
289
290           Bool_t hasDecayed = kTRUE;
291           if (idDec[iproc] == 0) Decay2Body(mother);
292           else if (idDec[iproc] == 1) DalitzDecay(mother); 
293           else DecayPiK(mother,hasDecayed); 
294           if (!hasDecayed) continue; 
295           Bool_t isMu0Acc = KinematicSelection(fMu[0],1); 
296           Bool_t isMu1Acc = KinematicSelection(fMu[1],1); 
297           Bool_t isMuFromPiKAcc = kTRUE; 
298
299           if (idDec[iproc] == 2) isMuFromPiKAcc = (mother->GetPdgCode()>0) ? isMu0Acc : isMu1Acc;
300           // mother 
301           if ((idDec[iproc]  < 2 && (isMu0Acc || isMu1Acc)) || 
302               (idDec[iproc] == 2 && isMuFromPiKAcc)) { 
303             pdgPushed[npartPushed] = mother->GetPdgCode(); 
304             pxPushed[npartPushed] = mother->Px(); 
305             pyPushed[npartPushed] = mother->Py(); 
306             pzPushed[npartPushed] = mother->Pz();
307             ePushed[npartPushed] = mother->Energy(); 
308             npartPushed++; 
309             if (isMu0Acc && (idDec[iproc] < 2 || mother->GetPdgCode() > 0)) { 
310               pdgPushed[npartPushed] = fMu[0]->GetPdgCode(); 
311               pxPushed[npartPushed] = fMu[0]->Px(); 
312               pyPushed[npartPushed] = fMu[0]->Py(); 
313               pzPushed[npartPushed] = fMu[0]->Pz();
314               ePushed[npartPushed] = fMu[0]->Energy(); 
315               npartPushed++; 
316               nmuons++; 
317             }
318             
319             if (isMu1Acc && (idDec[iproc] < 2 || mother->GetPdgCode() < 0)) { 
320               pdgPushed[npartPushed] = fMu[1]->GetPdgCode(); 
321               pxPushed[npartPushed] = fMu[1]->Px(); 
322               pyPushed[npartPushed] = fMu[1]->Py(); 
323               pzPushed[npartPushed] = fMu[1]->Pz();
324               ePushed[npartPushed] = fMu[1]->Energy(); 
325               npartPushed++; 
326               nmuons++; 
327             }
328           } 
329         } // end if BR
330       } // end loop on multiplicity 
331     }  // end loop on process 
332     fHMultMu->Fill(nmuons); 
333   } // keep on generating until at least a muon is created in the event
334   
335   Int_t ntmother = 0, ntchild =0; 
336   for (Int_t ipart = 0; ipart < npartPushed; ipart++) { 
337     if (TMath::Abs(pdgPushed[ipart]) != 13) { // particle is not a muon, hence it's a mother
338       PushTrack(0,-1,pdgPushed[ipart],
339                 pxPushed[ipart],pyPushed[ipart],pzPushed[ipart],ePushed[ipart],
340                 origin0[0],origin0[1],origin0[2],0.,
341                 polar[0],polar[1],polar[2],
342                 kPPrimary,ntmother,1,11);
343       KeepTrack(ntmother); 
344     }
345     else { 
346       PushTrack(1,ntmother,pdgPushed[ipart],
347                 pxPushed[ipart],pyPushed[ipart],pzPushed[ipart],ePushed[ipart],
348                 origin0[0],origin0[1],origin0[2],0.,
349                 polar[0],polar[1],polar[2],
350                 kPDecay,ntchild,1,1);
351       KeepTrack(ntchild); 
352     }
353   }
354   SetHighWaterMark(ntchild); 
355   AliGenEventHeader* header = new AliGenEventHeader("LMR");
356   header->SetPrimaryVertex(fVertex);
357   header->SetNProduced(fNprimaries);
358   AddHeader(header); 
359 }
360
361 //------------------------------------------------------------------
362
363 void AliGenMUONLMR::Decay2Body(TParticle *mother){ 
364   // performs decay in two muons of the low mass resonances
365   Double_t md1 = fMu[0]->GetMass(); 
366   Int_t pdg = mother->GetPdgCode(); 
367   Double_t mres =0; 
368   // if mother is a rho, extract the mass from its line shape
369   // otherwise consider the resonance mass 
370   if (pdg == 113) mres = fRhoLineShape->GetRandom(); 
371   else mres = mother->GetCalcMass(); 
372   //  while (mres < md1 + md2) mres =  fDsigmaDm[res]->GetRandom();
373   // energies and momenta in rest frame 
374   Double_t e1 = mres / 2.;
375   Double_t p1 = TMath::Sqrt((e1 + md1)*(e1 - md1)); 
376   // orientation in decaying particle rest frame
377   Double_t costheta = gRandom->Rndm() * 2 - 1;
378   Double_t sintheta = TMath::Sqrt((1. + costheta)*(1. - costheta));
379   Double_t phi      = 2. * TMath::Pi() * gRandom->Rndm(); 
380   Double_t px1      = p1 * sintheta * TMath::Cos(phi); 
381   Double_t py1      = p1 * sintheta * TMath::Sin(phi); 
382   Double_t pz1      = p1 * costheta; 
383
384   // boost muons into lab frame 
385
386   TLorentzVector vmother, v1, v2;
387   //  TLorentzVector boosted1, boosted2;   
388   vmother.SetPxPyPzE(mother->Px(),mother->Py(),mother->Pz(),mother->Energy());
389   v1.SetPxPyPzE(px1,py1,pz1,e1); 
390   v2.SetPxPyPzE(-px1,-py1,-pz1,e1); 
391
392   TVector3 betaParent = (1./vmother.E())*vmother.Vect(); // beta = p/E
393   v1.Boost(betaParent);
394   v2.Boost(betaParent);
395
396 //   TLorentzVector vtot = v1 + v2; 
397 //   printf ("mother: %g   %g    %g     %g\n",vmother.Px(), vmother.Py(), vmother.Pz(), vmother.E());
398 //   printf ("vtot  : %g   %g    %g     %g\n",vtot.Px(), vtot.Py(), vtot.Pz(), vtot.E());
399
400   fMu[0]->SetMomentum(v1.Px(),v1.Py(),v1.Pz(),v1.E());
401   fMu[1]->SetMomentum(v2.Px(),v2.Py(),v2.Pz(),v2.E());
402
403
404 //------------------------------------------------------------------
405
406 void AliGenMUONLMR::DecayPiK(TParticle *mother, Bool_t &hasDecayed){ 
407   // performs decays of pions and kaons
408   Double_t md1 = fMu[0]->GetMass(); 
409   // extract the mass from the resonance's line shape
410   Double_t mres = mother->GetMass(); 
411   // choose the pi/k sign, assuming 50% probabilities for both signs
412   Int_t sign = (gRandom->Rndm() > 0.5) ? 1 : -1;
413   mother->SetPdgCode(sign * TMath::Abs(mother->GetPdgCode())); 
414
415   // energies and momenta in rest frame 
416   Double_t e1 = (mres*mres + md1*md1)/(2*mres);
417   Double_t p1 = TMath::Sqrt((e1 + md1)*(e1 - md1)); 
418   // orientation in decaying particle rest frame
419   Double_t costheta = gRandom->Rndm() * 2 - 1;
420   Double_t sintheta = TMath::Sqrt((1. + costheta)*(1. - costheta));
421   Double_t phi      = 2. * TMath::Pi() * gRandom->Rndm(); 
422   Double_t px1      = p1 * sintheta * TMath::Cos(phi); 
423   Double_t py1      = p1 * sintheta * TMath::Sin(phi); 
424   Double_t pz1      = p1 * costheta;  
425
426   // boost muons into lab frame 
427   TLorentzVector vmother, v1;
428   vmother.SetPxPyPzE(mother->Px(),mother->Py(),mother->Pz(),mother->Energy());
429   v1.SetPxPyPzE(px1,py1,pz1,e1); 
430
431   TVector3 betaParent = (1./vmother.E())*vmother.Vect(); // beta = p/E
432   v1.Boost(betaParent);  
433   if (mother->GetPdgCode()>0) fMu[0]->SetMomentum(v1.Px(),v1.Py(),v1.Pz(),v1.E());
434   else fMu[1]->SetMomentum(v1.Px(),v1.Py(),v1.Pz(),v1.E());
435
436   Int_t idmother = -1; 
437   if (TMath::Abs(mother->GetPdgCode())== 211) idmother = 0; 
438   if (TMath::Abs(mother->GetPdgCode())== 321) idmother = 1; 
439   Double_t gammaRes = mother->Energy()/mres;
440   Double_t zResCM = fDecay[idmother]->GetRandom();
441   Double_t zResLab = gammaRes*zResCM;  
442   if(zResLab > 0.938) hasDecayed = 0; // 0.938: distance from IP to absorber + lambda_i
443   else hasDecayed = 1;
444
445
446
447 //-------------------------------------------------------------------
448
449 void AliGenMUONLMR::DalitzDecay(TParticle *mother){
450   //
451   // perform dalitz decays of eta, omega and etaprime 
452   //
453   //in the rest frame of the virtual photon:
454   Double_t mres = mother->GetCalcMass(); 
455   Double_t mumass  = fMu[0]->GetMass(); 
456   Double_t md3  = 0;  // unless differently specified, third particle is a photon 
457   if (mother->GetPdgCode() == 223) md3 = 0.134977; // if mother is an omega, third particle is a pi0
458   Int_t index = -1; 
459   if (mother->GetPdgCode() == 221) index = 0;  // eta
460   else if (mother->GetPdgCode() == 223) index = 1; // omega  
461   else if (mother->GetPdgCode() == 331) index = 2; // etaPrime  
462   Int_t flag = 0; 
463   Double_t xd=0, mvirt2=0; 
464   Double_t countIt = 0;
465   while (flag==0) {  
466     xd       = fDalitz[index]->GetRandom(); 
467     mvirt2   = xd * mres * mres;   // mass of virtual photon 
468     // check kinematics 
469     if (mres - md3 > TMath::Sqrt(mvirt2) && TMath::Sqrt(mvirt2)/2. > mumass) flag=1;
470     if (++countIt>1E11) {
471       mvirt2 =  mres * mres * 0.998; 
472       break;
473     }
474   }  
475
476   //
477   //        Generate muons in virtual photon rest frame. 
478   //        z axis is the virt. photon direction (before boost)  
479   //
480
481   Double_t e1 = TMath::Sqrt(mvirt2)/2.; // energy of mu1 in the virtual photon frame
482   Double_t psquare = (e1 + mumass)*(e1 - mumass); 
483   if (psquare<0) {
484     printf("Error in AliGenMUONLMR::DalitzDecay: sqrt of psquare = %f put to 0\n",psquare); 
485     psquare = 0;
486   }
487   Double_t p1 = TMath::Sqrt(psquare);
488   //theta angle between the pos. muon and the virtual photon 
489   Double_t costheta = fCosTheta->GetRandom();
490   if (costheta>1)  costheta = 1; 
491   if (costheta<-1) costheta = -1; 
492   Double_t sintheta = TMath::Sqrt((1. + costheta)*(1. - costheta));
493   Double_t phi      = 2 * TMath::Pi() * gRandom->Rndm();
494   Double_t sinphi   = TMath::Sin(phi);
495   Double_t cosphi   = TMath::Cos(phi);
496
497   // fill 4-vectors of leptons in the virtual photon frame
498
499   Double_t px1 = p1*sintheta*cosphi; 
500   Double_t py1 = p1*sintheta*sinphi; 
501   Double_t pz1 = p1*costheta; 
502   Double_t px2 = -p1*sintheta*cosphi; 
503   Double_t py2 = -p1*sintheta*sinphi; 
504   Double_t pz2 = -p1*costheta; 
505   Double_t e2  = e1; 
506
507   fMu[0]->SetMomentum(px1,py1,pz1,e1); 
508   fMu[1]->SetMomentum(px2,py2,pz2,e2); 
509
510   // calculate components of non-dilepton in CMS of parent resonance 
511
512   Double_t e3 = (mres * mres + md3 * md3 - mvirt2) / (2.*mres);
513   Double_t psquare3 = (e3 + md3)*(e3 - md3); 
514   if (psquare3<0) {
515     printf("Error in AliGenMUONLMR::DalitzDecay: sqrt of psquare3 = %f put to 0\n",psquare3); 
516     psquare3 = 0;
517   }
518   Double_t p3 = TMath::Sqrt(psquare3);
519   Double_t costheta2 = 2.* gRandom->Rndm() - 1.;   // angle between virtual photon and resonance
520   if (costheta2>1)  costheta2 = 1; 
521   if (costheta2<-1) costheta2 = -1; 
522   Double_t sintheta2 = TMath::Sqrt((1. + costheta2)*(1. - costheta2));
523   Double_t phi2      = 2 * TMath::Pi() * gRandom->Rndm();
524   Double_t sinphi2   = TMath::Sin(phi2);
525   Double_t cosphi2   = TMath::Cos(phi2);
526   Double_t px3 = p3*sintheta2*cosphi2; 
527   Double_t py3 = p3*sintheta2*sinphi2; 
528   Double_t pz3 = p3*costheta2; 
529   TLorentzVector v3(px3,py3,pz3,e3); 
530
531   sintheta2 = -sintheta2;
532   cosphi2   = -cosphi2;
533   sinphi2   = -sinphi2;
534
535   Double_t px1new = px1*costheta2*cosphi2 - py1*sinphi2 + pz1*sintheta2*cosphi2; 
536   Double_t py1new = px1*costheta2*sinphi2 + py1*cosphi2 + pz1*sintheta2*sinphi2; 
537   Double_t pz1new =-px1*sintheta2                       + pz1*costheta2; 
538   Double_t px2new = px2*costheta2*cosphi2 - py2*sinphi2 + pz2*sintheta2*cosphi2; 
539   Double_t py2new = px2*costheta2*sinphi2 + py2*cosphi2 + pz2*sintheta2*sinphi2; 
540   Double_t pz2new =-px2*sintheta2                       + pz2*costheta2; 
541
542   fMu[0]->SetMomentum(px1new,py1new,pz1new,e1); 
543   fMu[1]->SetMomentum(px2new,py2new,pz2new,e2); 
544
545   Double_t evirt = mres - e3; 
546   Double_t pxvirt = -px3;
547   Double_t pyvirt = -py3;
548   Double_t pzvirt = -pz3;
549   TLorentzVector vvirt(pxvirt,pyvirt,pzvirt,evirt); 
550   TVector3 betaVirt = (1./evirt) * vvirt.Vect(); // virtual photon beta in res frame
551
552   TLorentzVector v1(px1,py1,pz1,e1); 
553   TLorentzVector v2(px2,py2,pz2,e2);
554
555   // boost the muons in the frame where the resonance is at rest 
556
557   v1.Boost(betaVirt); 
558   v2.Boost(betaVirt); 
559
560   // boost muons and third particle in lab frame
561
562   TLorentzVector vmother(mother->Px(), mother->Py(), mother->Pz(), mother->Energy());  
563   TVector3 resBetaLab = (1./vmother.E())*vmother.Vect(); // eta beta in lab frame
564   v1.Boost(resBetaLab); 
565   v2.Boost(resBetaLab); 
566   v3.Boost(resBetaLab); 
567   vvirt.Boost(resBetaLab); 
568
569   fMu[0]->SetMomentum(v1.Px(),v1.Py(),v1.Pz(),v1.E());
570   fMu[1]->SetMomentum(v2.Px(),v2.Py(),v2.Pz(),v2.E());
571 //   part3->SetMomentum(v3.Px(),v3.Py(),v3.Pz(),v3.E());
572
573 //   TLorentzVector vtot = v1 + v2 + v3; 
574 //   TLorentzVector vdimu = v1 + v2; 
575 //   printf ("mother: %g   %g    %g     %g\n",vmother.Px(), vmother.Py(), vmother.Pz(), vmother.E());
576 //   printf ("vtot  : %g   %g    %g     %g\n",vtot.Px(), vtot.Py(), vtot.Pz(), vtot.E());
577 //   printf ("vvirt : %g   %g    %g     %g\n",vvirt.Px(), vvirt.Py(), vvirt.Pz(), vvirt.E());
578 //   printf ("vdimu : %g   %g    %g     %g\n",vdimu.Px(), vdimu.Py(), vdimu.Pz(), vdimu.E());
579
580 }
581
582 //------------------------------------------------------------------
583
584 Double_t AliGenMUONLMR::FormFactor(Double_t q2, Int_t decay){ 
585   //  Calculates the form factor for Dalitz decays A->B+l+l
586   //  Returns: |F(q^2)|^2
587   //
588   //  References: L.G. Landsberg, Physics Reports 128 No.6 (1985) 301-376. 
589  
590   Double_t ff2, mass2;
591   Double_t n2, n4, m2; 
592   // Lepton-G
593   
594   Double_t lambda2inv = 0;
595   switch (decay) { 
596   case 0:   // eta -> mu mu gamma  
597   // eta   -> l+ l- gamma: pole approximation
598     lambda2inv = 1.95; 
599     mass2 = fParticle[kEtaLMR]->GetMass() * fParticle[kEtaLMR]->GetMass(); 
600     if (q2 < mass2) ff2 = TMath::Power(1./(1.-lambda2inv*q2),2);
601     else ff2 = 0; 
602     break;
603   case 1:   // omega -> mu mu pi0 
604     // omega -> l+ l- pi0: pole approximation
605     mass2 = fParticle[kOmegaLMR]->GetMass() * fParticle[kOmegaLMR]->GetMass(); 
606     lambda2inv = 2.26; 
607     if (q2 < mass2) ff2 = TMath::Power(1./(1.-lambda2inv*q2),2);
608     else ff2 = 0; 
609     break;
610   case 2:   // etaPrime -> mu mu gamma 
611     mass2 = fParticle[kEtaPrimeLMR]->GetMass() * fParticle[kEtaPrimeLMR]->GetMass(); 
612     // eta'  -> l+ l- gamma: Breit-Wigner fitted to data
613     n2 = 0.764 * 0.764; 
614     n4 = n2 * n2; 
615     m2 = 0.1020 * 0.1020;
616     if (q2 < mass2) ff2 = n4 / (TMath::Power(n2-q2,2) + m2 * n2); 
617     else ff2 = 0; 
618     break;
619   default:
620     printf ("FormFactor: Decay not found\n"); 
621     return 0; 
622     break; 
623   }
624   return ff2; 
625 }
626
627 //____________________________________________________________
628
629 Double_t AliGenMUONLMR::RhoLineShapeNew(Double_t *x, Double_t* /*para*/){
630   //new parameterization implemented by Hiroyuki Sako (GSI)
631   Double_t mass = *x;
632   double r, GammaTot;
633   Double_t mRho    = TDatabasePDG::Instance()->GetParticle("rho0")->Mass();
634   Double_t mPi     = TDatabasePDG::Instance()->GetParticle("pi0")->Mass();
635   Double_t mMu     = TDatabasePDG::Instance()->GetParticle("mu-")->Mass();
636   Double_t Gamma0  = TDatabasePDG::Instance()->GetParticle("rho0")->Width();
637
638   const double Norm = 0.0744416*1.01;  
639   // 0.0744416 at m = 0.72297
640   // is the max number with Norm=1 (for rho)
641   
642   double mThreshold = 2.*mPi;
643
644   const double T = 0.170; // Assumption of pi+ temperature [GeV/c^2]
645   //const double T = 0.11; // Taken from fit to pi+ temperature [GeV/c^2]
646   // with Reference: LEBC-EHS collab., Z. Phys. C 50 (1991) 405
647
648   if (mass < mThreshold) {
649     r = 0.;
650     return r;
651   }
652
653   double k = sqrt(0.25*mass*mass-(mThreshold/2)*(mThreshold/2));
654   double k0 = sqrt(0.25*mRho*mRho-(mThreshold/2)*(mThreshold/2));
655
656   GammaTot = (k/k0)*(k/k0)*(k/k0)*(mRho/mass)*(mRho/mass)*Gamma0;
657
658   double FormFactor2 = 1/((mass*mass-mRho*mRho)*(mass*mass-mRho*mRho)+
659                           mass*mass*GammaTot*GammaTot);
660
661   r = pow(mass,1.5)*pow((1-mThreshold*mThreshold/(mass*mass)),1.5)*
662     ((mass*mass+2*mMu*mMu)/(mass*mass))*(pow((mass*mass-4*mMu*mMu),0.5)/mass)*FormFactor2
663     *exp(-mass/T)/Norm;
664
665   return r;
666 }
667
668