1 //#include <TClonesArray.h>
3 #include <TDatabasePDG.h>
6 #include "AliGenMUONLMR.h"
9 #include "AliGenEventHeader.h"
11 ClassImp(AliGenMUONLMR)
13 AliGenMUONLMR::AliGenMUONLMR () : AliGenMC(), fNMuMin(2), fGenSingleProc(-1) {
15 // default constructor
17 // initialize pt and y distributions according to a fit to
18 // Pythia simulation at sqrt(s) = 7 TeV
19 printf ("Using AliGenMUONLMR as generator\n");
20 for (Int_t ipart=0; ipart < fgkNpart; ipart++) fScaleMult[ipart] = 1;
21 fScaleMult[kPionLMR] = 0; // set pion multiplicity to zero
22 fScaleMult[kKaonLMR] = 0; // set kaon multiplicity to zero
23 Int_t pdg[7] = {211, 321, 221, 113, 223, 333, 331};
24 char* fptname[7] = {"fPtPion","fPtKaon","fPtEta","fPtRho","fPtOmega","fPtPhi","fPtEtaPrime"};
25 char* fyname[7] = {"fYPion","fYKaon","fYEta","fYRho","fYOmega","fYPhi","fYEtaPrime"};
26 char* fnname[7] = {"fMultPion","fMultKaon","fMultEta","fMultRho","fMultOmega","fMultPhi","fMultEtaPrime"};
27 char* fdname[2] = {"fDecPion","fDecKaon"};
28 Double_t ptparam[7][3] = { {1,0.427,2.52}, // pions from Pythia
29 {1,0.58,2.57}, // kaons from Pythia
30 {1,0.641,2.62}, // eta from Pythia
31 {1,1.2,2.5}, // rho+omega from ALICE muon
32 {1,1.2,2.5}, // rho+omega from ALICE muon
33 {1,1.03,2.5}, // phi from ALICE muon
34 {1,0.72,2.5}}; // etaPrime from Pythia
36 Double_t yparam[7][3] = { {1, 0.8251, 3.657}, // pions from pythia
37 {1, 1.83, 2.698}, // kaons from pythia
38 {1, 1.169, 3.282}, // eta from pythia
39 {1, 1.234, 3.264}, // rho from pythia
40 {1, 1.311, 3.223}, // omega from pythia
41 {1, 2.388, 2.129}, // phi from pythia
42 {1, 1.13,3.3}}; // eta prime from pythia
44 // multiplicity parameters from pythia
45 Double_t nparam[7][9] = { {353.582, 6.76263, 1.66979, 998.445, 9.73281, 12.6704, 175.187, 29.08, 40.2531},
46 {1.e4, 0.2841, 0,0,0,0,0,0,0},
47 {1.e4, 0.2647, 0,0,0,0,0,0,0},
48 {7055, 0.1786, 0,0,0,0,0,0,0},
49 {7500, 0.1896, 0,0,0,0,0,0,0},
50 {5.e04, 1.167, 0,0,0,0,0,0,0},
51 {2.9e04,0.714, 0,0,0,0,0,0,0}};
53 Double_t ctau[2] = {7.8045, 3.712};
55 for (Int_t i=0; i<fgkNpart; i++) {
58 fMult[i] = new TF1(fnname[i],"[0]*exp(-[1]*x)",0,30);
59 fMult[i]->SetParameters(nparam[i][0],nparam[i][1]);
62 fMult[i] = new TF1(fnname[i],"gaus(0)+gaus(3)+gaus(6)",0,150);
63 for (Int_t j=0; j<9; j++) fMult[i]->SetParameter(j,nparam[i][j]);
66 fPt[i] = new TF1(fptname[i],AliGenMUONLMR::PtDistr,0,20,3);
67 fPt[i]->SetParameters(ptparam[i][0], ptparam[i][1], ptparam[i][2]);
68 fY[i] = new TF1(fyname[i],AliGenMUONLMR::YDistr,-10,10,3);
69 fY[i]->SetParameters(yparam[i][0], yparam[i][1], yparam[i][2]);
72 for(Int_t i = 0; i<2; i++){
73 fDecay[i] = new TF1(fdname[i],"exp(-x/[0])",0,150);
74 fDecay[i]->SetParameter(0,ctau[i]);
77 for (Int_t ipart = 0; ipart < fgkNpart; ipart++) {
78 fParticle[ipart] = new TParticle();
79 fParticle[ipart]->SetPdgCode(fPDG[ipart]);
82 TDatabasePDG *pdgdb = TDatabasePDG::Instance();
83 Double_t mumass = pdgdb->GetParticle(13)->Mass();
84 fMu[0] = new TParticle();
85 fMu[0]->SetPdgCode(-13);
86 fMu[0]->SetCalcMass(mumass);
87 fMu[1] = new TParticle();
88 fMu[1]->SetPdgCode(13);
89 fMu[1]->SetCalcMass(mumass);
91 // function for polarized theta distributions
92 fCosTheta = new TF1 ("fCosTheta","1+[0]*x*x",-1,1);
93 fCosTheta->SetParameter(0,1);
97 Double_t xmin = 0, xmax = 2;
98 fDalitz[0] = new TH1F("hDalitzEta","",nbins,xmin,xmax);
99 fDalitz[1] = new TH1F("hDalitzOmega","",nbins,xmin,xmax);
100 fDalitz[2] = new TH1F("hDalitzEtaPrime","",nbins,xmin,xmax);
102 Double_t meta = pdgdb->GetParticle("eta")->Mass();
103 Double_t momega = pdgdb->GetParticle("omega")->Mass();
104 Double_t metaPrime = pdgdb->GetParticle("eta'")->Mass();
105 Double_t mpi0 = pdgdb->GetParticle("pi0")->Mass();
106 Double_t md3 = 0, mres = 0;
108 for (Int_t index = 0; index < 3; index++) {
113 else if (index == 1) {
117 else if (index == 2) {
121 Double_t delta = md3 * md3 / (mres * mres);
122 Double_t epsilon = mumass * mumass / (mres * mres);
123 Int_t nbins = fDalitz[index]->GetNbinsX();
124 Double_t xmin = fDalitz[index]->GetXaxis()->GetXmin();
125 Double_t deltax = fDalitz[index]->GetBinWidth(1);
126 Double_t xd = xmin - deltax/2.;
127 for (Int_t ibin = 0; ibin< nbins; ibin++) {
130 if (xd > 4. *epsilon) {
131 Double_t bracket = TMath::Power(1. + xd/(1. - delta),2)
132 - 4. * xd / ((1. - delta) * (1. - delta));
134 dalval = TMath::Power(bracket,1.5) /xd *
135 TMath::Sqrt(1 - 4 * epsilon / xd) * (1 + 2 * epsilon / xd) *
136 FormFactor(xd * mres * mres, index);
137 fDalitz[index]->Fill(xd,dalval);
143 fRhoLineShape = new TF1("fRhoLineShape",RhoLineShapeNew,0,2,2);
144 fHMultMu = new TH1D("fHMultMu","Muon multiplicity",20,-0.5,19.5);
145 fHNProc = new TH1D("fHNProc","Number of gen. evts. per process in 4 pi",9,-0.5,8.5);
148 //-----------------------------------------------------------
150 AliGenMUONLMR::~AliGenMUONLMR()
152 // Default destructor
153 for (Int_t i=0; i<7; i++) {
160 for (Int_t i=0; i<2; i++) {
165 for (Int_t i=0; i<3; i++) delete fDalitz[i];
167 delete fCosTheta; fCosTheta = 0;
168 delete fRhoLineShape; fRhoLineShape = 0;
169 delete fHMultMu; fHMultMu = 0;
170 delete fHNProc; fHNProc = 0;
173 //-----------------------------------------------------------
175 void AliGenMUONLMR::FinishRun(){
176 // save some histograms to an output file
177 Int_t nbins = fHNProc->GetNbinsX();
178 for (Int_t ibin=1; ibin <= nbins; ibin++) printf ("ibin = %d nEvProc = %g\n",
179 ibin,fHNProc->GetBinContent(ibin));
180 TFile *fout = new TFile("AliGenMUONLMR_histos.root","recreate");
186 //-----------------------------------------------------------
188 Double_t AliGenMUONLMR::YDistr(Double_t *px, Double_t *par){
189 // function for rapidity distribution: plateau at par[0] +
190 // gaussian tails centered at par[1] and with par[2]=sigma
191 Double_t x = TMath::Abs(px[0]);
193 if (x<par[1]) func = par[0];
195 Double_t z = (x-par[1])/(par[2]);
196 func = par[0] * TMath::Exp(-0.5 * z * z);
201 //-----------------------------------------------------------
203 Double_t AliGenMUONLMR::PtDistr(Double_t *px, Double_t *par){
204 // pt distribution: power law
206 Double_t func = par[0] * x / TMath::Power((1+(x/par[1])*(x/par[1])),par[2]);
210 //-----------------------------------------------------------
212 void AliGenMUONLMR::Generate() {
214 // generate the low mass resonances and their decays according to
215 // the multiplicity parameterized by pythia and BR from PDG
216 // rapidity distributions parametrized from pythia
217 // pt distributions from data (or pythia for etaprime)
219 Double_t pxPushed[100], pyPushed[100], pzPushed[100], ePushed[100];
220 Int_t nmuons = -1, npartPushed = 0, pdgPushed[100];
221 Double_t polar[3]= {0,0,0}; // Polarisation of the parent particle (for GEANT tracking)
222 Double_t origin0[3]; // Origin of the generated parent particle (for GEANT tracking)
223 // Calculating vertex position per event
224 for (Int_t j=0;j<3;j++) origin0[j]=fOrigin[j];
225 if(fVertexSmear==kPerEvent) {
227 for (Int_t j=0;j<3;j++) origin0[j]=fVertex[j];
230 printf ("In Generate()\n");
232 TDatabasePDG* pdg = TDatabasePDG::Instance();
234 Double_t pt, y, phi, mass, px, py, pz, ene, mt;
236 const Int_t nproc = 9;
237 Int_t idRes[nproc] = {kEtaLMR, kEtaLMR, kRhoLMR, kOmegaLMR, kOmegaLMR, kPhiLMR, kEtaPrimeLMR, kPionLMR, kKaonLMR};
238 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};
239 // Double_t BR[nproc] = {1, 1, 1, 1, 1, 1, 1, 1, 1};
240 Int_t idDec[nproc] = {0, 1, 0, 0, 1, 0, 1, 2, 2}; // 0:2body, 1:Dalitz, 2:pi/K
241 Int_t mult[nproc] = {0,0,0,0,0,0,0,0,0};
243 while (nmuons < fNMuMin) {
247 for (Int_t iproc=0; iproc<nproc; iproc++) {
248 if (fGenSingleProc == -1) {
249 mult[iproc] = Int_t(fMult[idRes[iproc]]->GetRandom()*fScaleMult[idRes[iproc]]);
252 if (iproc==fGenSingleProc) {
263 if (fGenSingleProc == -1) {
268 for (Int_t iproc = 0; iproc < nproc; iproc++) {
269 // printf ("Multiplicity for process %d is %d\n",iproc,mult[iproc]);
270 for (Int_t imult=0; imult<mult[iproc]; imult++) {
271 if (gRandom->Rndm() < BR[iproc]) {
272 fHNProc->Fill(iproc);
273 Int_t ipart = idRes[iproc];
274 pt = fPt[ipart]->GetRandom();
275 y = fY[ipart]->GetRandom();
276 phi = gRandom->Rndm() * 2 * TMath::Pi();
277 mass = pdg->GetParticle(fPDG[ipart])->Mass();
278 px = pt * TMath::Cos(phi);
279 py = pt * TMath::Sin(phi);
280 mt = TMath::Sqrt(pt * pt + mass * mass);
281 pz = mt * TMath::SinH(y);
282 ene = mt * TMath::CosH(y);
284 mother = fParticle[ipart];
285 mother->SetMomentum(px,py,pz,ene);
286 mother->SetCalcMass(mass);
287 if (!KinematicSelection(mother,0)) continue;
289 Bool_t hasDecayed = kTRUE;
290 if (idDec[iproc] == 0) Decay2Body(mother);
291 else if (idDec[iproc] == 1) DalitzDecay(mother);
292 else DecayPiK(mother,hasDecayed);
293 if (!hasDecayed) continue;
294 Bool_t isMu0Acc = KinematicSelection(fMu[0],1);
295 Bool_t isMu1Acc = KinematicSelection(fMu[1],1);
296 Bool_t isMuFromPiKAcc = kTRUE;
298 if (idDec[iproc] == 2) isMuFromPiKAcc = (mother->GetPdgCode()>0) ? isMu0Acc : isMu1Acc;
300 if ((idDec[iproc] < 2 && (isMu0Acc || isMu1Acc)) ||
301 (idDec[iproc] == 2 && isMuFromPiKAcc)) {
302 pdgPushed[npartPushed] = mother->GetPdgCode();
303 pxPushed[npartPushed] = mother->Px();
304 pyPushed[npartPushed] = mother->Py();
305 pzPushed[npartPushed] = mother->Pz();
306 ePushed[npartPushed] = mother->Energy();
308 if (isMu0Acc && (idDec[iproc] < 2 || mother->GetPdgCode() > 0)) {
309 pdgPushed[npartPushed] = fMu[0]->GetPdgCode();
310 pxPushed[npartPushed] = fMu[0]->Px();
311 pyPushed[npartPushed] = fMu[0]->Py();
312 pzPushed[npartPushed] = fMu[0]->Pz();
313 ePushed[npartPushed] = fMu[0]->Energy();
318 if (isMu1Acc && (idDec[iproc] < 2 || mother->GetPdgCode() < 0)) {
319 pdgPushed[npartPushed] = fMu[1]->GetPdgCode();
320 pxPushed[npartPushed] = fMu[1]->Px();
321 pyPushed[npartPushed] = fMu[1]->Py();
322 pzPushed[npartPushed] = fMu[1]->Pz();
323 ePushed[npartPushed] = fMu[1]->Energy();
329 } // end loop on multiplicity
330 } // end loop on process
331 fHMultMu->Fill(nmuons);
332 } // keep on generating until at least a muon is created in the event
334 Int_t ntmother = 0, ntchild =0;
335 for (Int_t ipart = 0; ipart < npartPushed; ipart++) {
336 if (TMath::Abs(pdgPushed[ipart]) != 13) { // particle is not a muon, hence it's a mother
337 PushTrack(0,-1,pdgPushed[ipart],
338 pxPushed[ipart],pyPushed[ipart],pzPushed[ipart],ePushed[ipart],
339 origin0[0],origin0[1],origin0[2],0.,
340 polar[0],polar[1],polar[2],
341 kPPrimary,ntmother,1,11);
345 PushTrack(1,ntmother,pdgPushed[ipart],
346 pxPushed[ipart],pyPushed[ipart],pzPushed[ipart],ePushed[ipart],
347 origin0[0],origin0[1],origin0[2],0.,
348 polar[0],polar[1],polar[2],
349 kPDecay,ntchild,1,1);
353 SetHighWaterMark(ntchild);
354 AliGenEventHeader* header = new AliGenEventHeader("LMR");
355 header->SetPrimaryVertex(fVertex);
356 header->SetNProduced(fNprimaries);
360 //------------------------------------------------------------------
362 void AliGenMUONLMR::Decay2Body(TParticle *mother){
363 // performs decay in two muons of the low mass resonances
364 Double_t md1 = fMu[0]->GetMass();
365 Int_t pdg = mother->GetPdgCode();
367 // if mother is a rho, extract the mass from its line shape
368 // otherwise consider the resonance mass
369 if (pdg == 113) mres = fRhoLineShape->GetRandom();
370 else mres = mother->GetCalcMass();
371 // while (mres < md1 + md2) mres = fDsigmaDm[res]->GetRandom();
372 // energies and momenta in rest frame
373 Double_t e1 = mres / 2.;
374 Double_t p1 = TMath::Sqrt((e1 + md1)*(e1 - md1));
375 // orientation in decaying particle rest frame
376 Double_t costheta = gRandom->Rndm() * 2 - 1;
377 Double_t sintheta = TMath::Sqrt((1. + costheta)*(1. - costheta));
378 Double_t phi = 2. * TMath::Pi() * gRandom->Rndm();
379 Double_t px1 = p1 * sintheta * TMath::Cos(phi);
380 Double_t py1 = p1 * sintheta * TMath::Sin(phi);
381 Double_t pz1 = p1 * costheta;
383 // boost muons into lab frame
385 TLorentzVector vmother, v1, v2;
386 // TLorentzVector boosted1, boosted2;
387 vmother.SetPxPyPzE(mother->Px(),mother->Py(),mother->Pz(),mother->Energy());
388 v1.SetPxPyPzE(px1,py1,pz1,e1);
389 v2.SetPxPyPzE(-px1,-py1,-pz1,e1);
391 TVector3 betaParent = (1./vmother.E())*vmother.Vect(); // beta = p/E
392 v1.Boost(betaParent);
393 v2.Boost(betaParent);
395 // TLorentzVector vtot = v1 + v2;
396 // printf ("mother: %g %g %g %g\n",vmother.Px(), vmother.Py(), vmother.Pz(), vmother.E());
397 // printf ("vtot : %g %g %g %g\n",vtot.Px(), vtot.Py(), vtot.Pz(), vtot.E());
399 fMu[0]->SetMomentum(v1.Px(),v1.Py(),v1.Pz(),v1.E());
400 fMu[1]->SetMomentum(v2.Px(),v2.Py(),v2.Pz(),v2.E());
403 //------------------------------------------------------------------
405 void AliGenMUONLMR::DecayPiK(TParticle *mother, Bool_t &hasDecayed){
406 // performs decays of pions and kaons
407 Double_t md1 = fMu[0]->GetMass();
408 // extract the mass from the resonance's line shape
409 Double_t mres = mother->GetMass();
410 // choose the pi/k sign, assuming 50% probabilities for both signs
411 Int_t sign = (gRandom->Rndm() > 0.5) ? 1 : -1;
412 mother->SetPdgCode(sign * TMath::Abs(mother->GetPdgCode()));
414 // energies and momenta in rest frame
415 Double_t e1 = (mres*mres + md1*md1)/(2*mres);
416 Double_t p1 = TMath::Sqrt((e1 + md1)*(e1 - md1));
417 // orientation in decaying particle rest frame
418 Double_t costheta = gRandom->Rndm() * 2 - 1;
419 Double_t sintheta = TMath::Sqrt((1. + costheta)*(1. - costheta));
420 Double_t phi = 2. * TMath::Pi() * gRandom->Rndm();
421 Double_t px1 = p1 * sintheta * TMath::Cos(phi);
422 Double_t py1 = p1 * sintheta * TMath::Sin(phi);
423 Double_t pz1 = p1 * costheta;
425 // boost muons into lab frame
426 TLorentzVector vmother, v1;
427 vmother.SetPxPyPzE(mother->Px(),mother->Py(),mother->Pz(),mother->Energy());
428 v1.SetPxPyPzE(px1,py1,pz1,e1);
430 TVector3 betaParent = (1./vmother.E())*vmother.Vect(); // beta = p/E
431 v1.Boost(betaParent);
432 if (mother->GetPdgCode()>0) fMu[0]->SetMomentum(v1.Px(),v1.Py(),v1.Pz(),v1.E());
433 else fMu[1]->SetMomentum(v1.Px(),v1.Py(),v1.Pz(),v1.E());
436 if (TMath::Abs(mother->GetPdgCode())== 211) idmother = 0;
437 if (TMath::Abs(mother->GetPdgCode())== 321) idmother = 1;
438 Double_t gammaRes = mother->Energy()/mres;
439 Double_t zResCM = fDecay[idmother]->GetRandom();
440 Double_t zResLab = gammaRes*zResCM;
441 if(zResLab > 0.938) hasDecayed = 0; // 0.938: distance from IP to absorber + lambda_i
446 //-------------------------------------------------------------------
448 void AliGenMUONLMR::DalitzDecay(TParticle *mother){
450 // perform dalitz decays of eta, omega and etaprime
452 //in the rest frame of the virtual photon:
453 Double_t mres = mother->GetCalcMass();
454 Double_t mumass = fMu[0]->GetMass();
455 Double_t md3 = 0; // unless differently specified, third particle is a photon
456 if (mother->GetPdgCode() == 223) md3 = 0.134977; // if mother is an omega, third particle is a pi0
458 if (mother->GetPdgCode() == 221) index = 0; // eta
459 else if (mother->GetPdgCode() == 223) index = 1; // omega
460 else if (mother->GetPdgCode() == 331) index = 2; // etaPrime
462 Double_t xd=0, mvirt2=0;
463 Double_t countIt = 0;
465 xd = fDalitz[index]->GetRandom();
466 mvirt2 = xd * mres * mres; // mass of virtual photon
468 if (mres - md3 > TMath::Sqrt(mvirt2) && TMath::Sqrt(mvirt2)/2. > mumass) flag=1;
469 if (++countIt>1E11) {
470 mvirt2 = mres * mres * 0.998;
476 // Generate muons in virtual photon rest frame.
477 // z axis is the virt. photon direction (before boost)
480 Double_t e1 = TMath::Sqrt(mvirt2)/2.; // energy of mu1 in the virtual photon frame
481 Double_t psquare = (e1 + mumass)*(e1 - mumass);
483 printf("Error in AliGenMUONLMR::DalitzDecay: sqrt of psquare = %f put to 0\n",psquare);
486 Double_t p1 = TMath::Sqrt(psquare);
487 //theta angle between the pos. muon and the virtual photon
488 Double_t costheta = fCosTheta->GetRandom();
489 if (costheta>1) costheta = 1;
490 if (costheta<-1) costheta = -1;
491 Double_t sintheta = TMath::Sqrt((1. + costheta)*(1. - costheta));
492 Double_t phi = 2 * TMath::Pi() * gRandom->Rndm();
493 Double_t sinphi = TMath::Sin(phi);
494 Double_t cosphi = TMath::Cos(phi);
496 // fill 4-vectors of leptons in the virtual photon frame
498 Double_t px1 = p1*sintheta*cosphi;
499 Double_t py1 = p1*sintheta*sinphi;
500 Double_t pz1 = p1*costheta;
501 Double_t px2 = -p1*sintheta*cosphi;
502 Double_t py2 = -p1*sintheta*sinphi;
503 Double_t pz2 = -p1*costheta;
506 fMu[0]->SetMomentum(px1,py1,pz1,e1);
507 fMu[1]->SetMomentum(px2,py2,pz2,e2);
509 // calculate components of non-dilepton in CMS of parent resonance
511 Double_t e3 = (mres * mres + md3 * md3 - mvirt2) / (2.*mres);
512 Double_t psquare3 = (e3 + md3)*(e3 - md3);
514 printf("Error in AliGenMUONLMR::DalitzDecay: sqrt of psquare3 = %f put to 0\n",psquare3);
517 Double_t p3 = TMath::Sqrt(psquare3);
518 Double_t costheta2 = 2.* gRandom->Rndm() - 1.; // angle between virtual photon and resonance
519 if (costheta2>1) costheta2 = 1;
520 if (costheta2<-1) costheta2 = -1;
521 Double_t sintheta2 = TMath::Sqrt((1. + costheta2)*(1. - costheta2));
522 Double_t phi2 = 2 * TMath::Pi() * gRandom->Rndm();
523 Double_t sinphi2 = TMath::Sin(phi2);
524 Double_t cosphi2 = TMath::Cos(phi2);
525 Double_t px3 = p3*sintheta2*cosphi2;
526 Double_t py3 = p3*sintheta2*sinphi2;
527 Double_t pz3 = p3*costheta2;
528 TLorentzVector v3(px3,py3,pz3,e3);
530 sintheta2 = -sintheta2;
534 Double_t px1new = px1*costheta2*cosphi2 - py1*sinphi2 + pz1*sintheta2*cosphi2;
535 Double_t py1new = px1*costheta2*sinphi2 + py1*cosphi2 + pz1*sintheta2*sinphi2;
536 Double_t pz1new =-px1*sintheta2 + pz1*costheta2;
537 Double_t px2new = px2*costheta2*cosphi2 - py2*sinphi2 + pz2*sintheta2*cosphi2;
538 Double_t py2new = px2*costheta2*sinphi2 + py2*cosphi2 + pz2*sintheta2*sinphi2;
539 Double_t pz2new =-px2*sintheta2 + pz2*costheta2;
541 fMu[0]->SetMomentum(px1new,py1new,pz1new,e1);
542 fMu[1]->SetMomentum(px2new,py2new,pz2new,e2);
544 Double_t evirt = mres - e3;
545 Double_t pxvirt = -px3;
546 Double_t pyvirt = -py3;
547 Double_t pzvirt = -pz3;
548 TLorentzVector vvirt(pxvirt,pyvirt,pzvirt,evirt);
549 TVector3 betaVirt = (1./evirt) * vvirt.Vect(); // virtual photon beta in res frame
551 TLorentzVector v1(px1,py1,pz1,e1);
552 TLorentzVector v2(px2,py2,pz2,e2);
554 // boost the muons in the frame where the resonance is at rest
559 // boost muons and third particle in lab frame
561 TLorentzVector vmother(mother->Px(), mother->Py(), mother->Pz(), mother->Energy());
562 TVector3 resBetaLab = (1./vmother.E())*vmother.Vect(); // eta beta in lab frame
563 v1.Boost(resBetaLab);
564 v2.Boost(resBetaLab);
565 v3.Boost(resBetaLab);
566 vvirt.Boost(resBetaLab);
568 fMu[0]->SetMomentum(v1.Px(),v1.Py(),v1.Pz(),v1.E());
569 fMu[1]->SetMomentum(v2.Px(),v2.Py(),v2.Pz(),v2.E());
570 // part3->SetMomentum(v3.Px(),v3.Py(),v3.Pz(),v3.E());
572 // TLorentzVector vtot = v1 + v2 + v3;
573 // TLorentzVector vdimu = v1 + v2;
574 // printf ("mother: %g %g %g %g\n",vmother.Px(), vmother.Py(), vmother.Pz(), vmother.E());
575 // printf ("vtot : %g %g %g %g\n",vtot.Px(), vtot.Py(), vtot.Pz(), vtot.E());
576 // printf ("vvirt : %g %g %g %g\n",vvirt.Px(), vvirt.Py(), vvirt.Pz(), vvirt.E());
577 // printf ("vdimu : %g %g %g %g\n",vdimu.Px(), vdimu.Py(), vdimu.Pz(), vdimu.E());
581 //------------------------------------------------------------------
583 Double_t AliGenMUONLMR::FormFactor(Double_t q2, Int_t decay){
584 // Calculates the form factor for Dalitz decays A->B+l+l
585 // Returns: |F(q^2)|^2
587 // References: L.G. Landsberg, Physics Reports 128 No.6 (1985) 301-376.
593 Double_t lambda2inv = 0;
595 case 0: // eta -> mu mu gamma
596 // eta -> l+ l- gamma: pole approximation
598 mass2 = fParticle[kEtaLMR]->GetMass() * fParticle[kEtaLMR]->GetMass();
599 if (q2 < mass2) ff2 = TMath::Power(1./(1.-lambda2inv*q2),2);
602 case 1: // omega -> mu mu pi0
603 // omega -> l+ l- pi0: pole approximation
604 mass2 = fParticle[kOmegaLMR]->GetMass() * fParticle[kOmegaLMR]->GetMass();
606 if (q2 < mass2) ff2 = TMath::Power(1./(1.-lambda2inv*q2),2);
609 case 2: // etaPrime -> mu mu gamma
610 mass2 = fParticle[kEtaPrimeLMR]->GetMass() * fParticle[kEtaPrimeLMR]->GetMass();
611 // eta' -> l+ l- gamma: Breit-Wigner fitted to data
614 m2 = 0.1020 * 0.1020;
615 if (q2 < mass2) ff2 = n4 / (TMath::Power(n2-q2,2) + m2 * n2);
619 printf ("FormFactor: Decay not found\n");
626 //____________________________________________________________
628 Double_t AliGenMUONLMR::RhoLineShapeNew(Double_t *x, Double_t *para){
629 //new parameterization implemented by Hiroyuki Sako (GSI)
632 Double_t mRho = TDatabasePDG::Instance()->GetParticle("rho0")->Mass();
633 Double_t mPi = TDatabasePDG::Instance()->GetParticle("pi0")->Mass();
634 Double_t mMu = TDatabasePDG::Instance()->GetParticle("mu-")->Mass();
635 Double_t Gamma0 = TDatabasePDG::Instance()->GetParticle("rho0")->Width();
637 const double Norm = 0.0744416*1.01;
638 // 0.0744416 at m = 0.72297
639 // is the max number with Norm=1 (for rho)
641 double mThreshold = 2.*mPi;
643 const double T = 0.170; // Assumption of pi+ temperature [GeV/c^2]
644 //const double T = 0.11; // Taken from fit to pi+ temperature [GeV/c^2]
645 // with Reference: LEBC-EHS collab., Z. Phys. C 50 (1991) 405
647 if (mass < mThreshold) {
652 double k = sqrt(0.25*mass*mass-(mThreshold/2)*(mThreshold/2));
653 double k0 = sqrt(0.25*mRho*mRho-(mThreshold/2)*(mThreshold/2));
655 GammaTot = (k/k0)*(k/k0)*(k/k0)*(mRho/mass)*(mRho/mass)*Gamma0;
657 double FormFactor2 = 1/((mass*mass-mRho*mRho)*(mass*mass-mRho*mRho)+
658 mass*mass*GammaTot*GammaTot);
660 r = pow(mass,1.5)*pow((1-mThreshold*mThreshold/(mass*mass)),1.5)*
661 ((mass*mass+2*mMu*mMu)/(mass*mass))*(pow((mass*mass-4*mMu*mMu),0.5)/mass)*FormFactor2