1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 // Class to generate particles from using paramtrized pT and y distributions.
19 // Distributions are obtained from pointer to object of type AliGenLib.
20 // (For example AliGenMUONlib)
21 // Decays are performed using Pythia.
22 // andreas.morsch@cern.ch
25 #include <TClonesArray.h>
26 #include <TDatabasePDG.h>
29 #include <TLorentzVector.h>
31 #include <TParticle.h>
32 #include <TParticlePDG.h>
34 #include <TVirtualMC.h>
36 #include "AliDecayer.h"
37 #include "AliGenMUONlib.h"
38 #include "AliGenParam.h"
41 #include "AliGenEventHeader.h"
45 //------------------------------------------------------------
49 <img src="picts/AliGenParam.gif">
53 //____________________________________________________________
54 AliGenParam::AliGenParam()
74 // Default constructor
76 //____________________________________________________________
77 AliGenParam::AliGenParam(Int_t npart, const AliGenLib * Library, Int_t param, const char* tname)
79 fPtParaFunc(Library->GetPt(param, tname)),
80 fYParaFunc (Library->GetY (param, tname)),
81 fIpParaFunc(Library->GetIp(param, tname)),
82 fV2ParaFunc(Library->GetV2(param, tname)),
98 // Constructor using number of particles parameterisation id and library
100 fTitle= "Particle Generator using pT and y parameterisation";
104 //____________________________________________________________
105 AliGenParam::AliGenParam(Int_t npart, Int_t param, const char* tname, const char* name):
126 // Constructor using parameterisation id and number of particles
129 fTitle= "Particle Generator using pT and y parameterisation";
131 AliGenLib* pLibrary = new AliGenMUONlib();
132 fPtParaFunc = pLibrary->GetPt(param, tname);
133 fYParaFunc = pLibrary->GetY (param, tname);
134 fIpParaFunc = pLibrary->GetIp(param, tname);
135 fV2ParaFunc = pLibrary->GetV2(param, tname);
139 for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
142 SetChildMomentumRange();
145 SetChildThetaRange();
147 //____________________________________________________________
149 AliGenParam::AliGenParam(Int_t npart, Int_t param,
150 Double_t (*PtPara) (const Double_t*, const Double_t*),
151 Double_t (*YPara ) (const Double_t* ,const Double_t*),
152 Double_t (*V2Para) (const Double_t* ,const Double_t*),
153 Int_t (*IpPara) (TRandom *))
176 // Gines Martinez 1/10/99
178 fTitle= "Particle Generator using pT and y parameterisation";
182 for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
185 SetChildMomentumRange();
188 SetChildThetaRange();
191 //____________________________________________________________
192 AliGenParam::~AliGenParam()
201 //-------------------------------------------------------------------
202 TVector3 AliGenParam::OrthogonalVector(TVector3 &inVec){
203 double abc[]={inVec.x(), inVec.y(), inVec.z()};
204 double xyz[]={1,1,1};
207 for(int i=0; i<3; i++)
212 xyz[solvDim]=(-abc[(1+solvDim)%3]-abc[(2+solvDim)%3])/abc[(0+solvDim)%3];
214 TVector3 res(xyz[0],xyz[1],xyz[2]);
218 double AliGenParam::ScreenFunc1(double d){
220 return 21.12-4.184*log(d+0.952);
222 return 20.867-3.242*d+0.652*d*d;
225 double AliGenParam::ScreenFunc2(double d){
227 return 21.12-4.184*log(d+0.952);
229 return 20.209-1.93*d-0.086*d*d;
232 double AliGenParam::EnergyFraction(double Z, double E){
233 double e0=0.000511/E;
235 double dmin=ScreenVar(Z,e0,0.5);
236 double fcZ=(aZ*aZ)*(1/(1+aZ*aZ)+0.20206-0.0368*aZ*aZ+0.0083*aZ*aZ*aZ);
237 double Fz=8*log(Z)/3;
240 double dmax=exp((42.24-Fz)/8.368)-0.952;
242 return fRandom->Uniform(e0,0.5);
244 double e1=0.5-0.5*sqrt(1-dmin/dmax);
245 double emin=TMath::Max(e0,e1);
247 double normval=1/(0.5*(ScreenFunc1(dmin)-0.5*Fz)+0.1666667*(ScreenFunc2(dmin)-0.5*Fz));
249 double y=fRandom->Uniform(emin,1);
250 double eps=y*(0.38453+y*(0.10234+y*(0.026072+y*(0.014367-0.027313*y)))); //inverse to the enveloping cumulative probability distribution
251 double val=fRandom->Uniform(0,2.12*(eps*eps-eps)+1.53); //enveloping probability density
252 double d=ScreenVar(Z,e0,eps);
253 double bh=((eps*eps)+(1-eps)*(1-eps))*(ScreenFunc1(d)-0.5*Fz)+0.6666667*eps*(1-eps)*(ScreenFunc2(d)-0.5*Fz);
260 double AliGenParam::PolarAngle(double E){
264 double u=-8*log(rand[1]*rand[2])/5;
265 if(!(rand[0]<9.0/36))
270 Int_t AliGenParam::ForceGammaConversion(TClonesArray *particles, Int_t nPart)
272 //based on: http://geant4.cern.ch/G4UsersDocuments/UsersGuides/PhysicsReferenceManual/html/node27.html
273 // and: http://geant4.cern.ch/G4UsersDocuments/UsersGuides/PhysicsReferenceManual/html/node58.html
274 Int_t nPartNew=nPart;
275 for(int iPart=0; iPart<nPart; iPart++){
276 TParticle *gamma = (TParticle *) particles->At(iPart);
277 if(gamma->GetPdgCode()!=22) continue;
279 TVector3 gammaV3(gamma->Px(),gamma->Py(),gamma->Pz());
280 Float_t az=fRandom->Uniform(TMath::Pi()*2);
281 double frac=EnergyFraction(1,gamma->Energy());
282 double Ee1=frac*gamma->Energy();
283 double Ee2=(1-frac)*gamma->Energy();
284 double Pe1=Ee1;//sqrt(Ee1*Ee1-0.000511*0.000511);
285 double Pe2=Ee2;//sqrt(Ee2*Ee2-0.000511*0.000511);
287 TVector3 rotAxis(OrthogonalVector(gammaV3));
288 rotAxis.Rotate(az,gammaV3);
289 TVector3 e1V3(gammaV3);
290 e1V3.Rotate(PolarAngle(Ee1),rotAxis);
293 TVector3 e2V3(gammaV3);
294 e2V3.Rotate(-PolarAngle(Ee2),rotAxis);
297 // gamma = new TParticle(*gamma);
298 // particles->RemoveAt(iPart);
299 gamma->SetFirstDaughter(nPartNew+1);
300 gamma->SetLastDaughter(nPartNew+2);
301 // new((*particles)[iPart]) TParticle(*gamma);
305 gamma->ProductionVertex(vtx);
306 new((*particles)[nPartNew]) TParticle(11, gamma->GetStatusCode(), iPart+1, -1, 0, 0, TLorentzVector(e1V3,Ee1), vtx);
308 new((*particles)[nPartNew]) TParticle(-11, gamma->GetStatusCode(), iPart+1, -1, 0, 0, TLorentzVector(e2V3,Ee2), vtx);
311 // particles->Compress();
312 return particles->GetEntriesFast();
315 //____________________________________________________________
316 void AliGenParam::Init()
320 if (gMC) fDecayer = gMC->GetDecayer();
323 <img src="picts/AliGenParam.gif">
327 snprintf(name, 256, "pt-parameterisation for %s", GetName());
329 if (fPtPara) fPtPara->Delete();
330 fPtPara = new TF1(name, fPtParaFunc, fPtMin, fPtMax,0);
331 gROOT->GetListOfFunctions()->Remove(fPtPara);
332 // Set representation precision to 10 MeV
333 Int_t npx= Int_t((fPtMax - fPtMin) / fDeltaPt);
335 fPtPara->SetNpx(npx);
337 snprintf(name, 256, "y-parameterisation for %s", GetName());
338 if (fYPara) fYPara->Delete();
339 fYPara = new TF1(name, fYParaFunc, fYMin, fYMax, 0);
340 gROOT->GetListOfFunctions()->Remove(fYPara);
342 snprintf(name, 256, "v2-parameterisation for %s", GetName());
343 if (fV2Para) fV2Para->Delete();
344 fV2Para = new TF1(name, fV2ParaFunc, fPtMin, fPtMax, 0);
345 // fV2Para = new TF1(name, "2*[0]/(1+TMath::Exp([1]*([2]-x)))-[0]", fPtMin, fPtMax);
346 // fV2Para->SetParameter(0, 0.236910);
347 // fV2Para->SetParameter(1, 1.71122);
348 // fV2Para->SetParameter(2, 0.0827617);
349 //gROOT->GetListOfFunctions()->Remove(fV2Para); //TR: necessary?
351 snprintf(name, 256, "dNdPhi for %s", GetName());
352 if (fdNdPhi) fdNdPhi->Delete();
353 fdNdPhi = new TF1(name, "1+2*[0]*TMath::Cos(2*(x-[1]))", fPhiMin, fPhiMax);
354 //gROOT->GetListOfFunctions()->Remove(fdNdPhi); //TR: necessary?
356 snprintf(name, 256, "pt-for-%s", GetName());
357 TF1 ptPara(name ,fPtParaFunc, 0, 15, 0);
358 snprintf(name, 256, "y-for-%s", GetName());
359 TF1 yPara(name, fYParaFunc, -6, 6, 0);
366 fdNdy0=fYParaFunc(&y1,&y2);
368 // Integral over generation region
369 #if ROOT_VERSION_CODE < ROOT_VERSION(5,99,0)
370 Float_t intYS = yPara.Integral(fYMin, fYMax,(Double_t*) 0x0,1.e-6);
371 Float_t intPt0 = ptPara.Integral(0,15,(Double_t *) 0x0,1.e-6);
372 Float_t intPtS = ptPara.Integral(fPtMin,fPtMax,(Double_t*) 0x0,1.e-6);
374 Float_t intYS = yPara.Integral(fYMin, fYMax,1.e-6);
375 Float_t intPt0 = ptPara.Integral(0,15,1.e-6);
376 Float_t intPtS = ptPara.Integral(fPtMin,fPtMax,1.e-6);
378 Float_t phiWgt=(fPhiMax-fPhiMin)/2./TMath::Pi(); //TR: should probably be done differently in case of anisotropic phi...
379 if (fAnalog == kAnalog) {
380 fYWgt = intYS/fdNdy0;
381 fPtWgt = intPtS/intPt0;
382 fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
384 fYWgt = intYS/fdNdy0;
385 fPtWgt = (fPtMax-fPtMin)/intPt0;
386 fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
389 // particle decay related initialization
390 fDecayer->SetForceDecay(fForceDecay);
397 //____________________________________________________________
398 void AliGenParam::Generate()
401 // Generate 1 event (see Generate(Int_t ntimes) for details
405 //____________________________________________________________
406 void AliGenParam::GenerateN(Int_t ntimes)
409 // Generate ntimes*'npart' light and heavy mesons (J/Psi, upsilon or phi, Pion,
410 // Kaons, Etas, Omegas) and Baryons (proton, antiprotons, neutrons and
411 // antineutrons in the the desired theta, phi and momentum windows;
412 // Gaussian smearing on the vertex is done if selected.
413 // The decay of heavy mesons is done using lujet,
414 // and the childern particle are tracked by GEANT
415 // However, light mesons are directly tracked by GEANT
416 // setting fForceDecay = nodecay (SetForceDecay(nodecay))
419 // Reinitialize decayer
420 fDecayer->SetForceDecay(fForceDecay);
424 Float_t polar[3]= {0,0,0}; // Polarisation of the parent particle (for GEANT tracking)
425 Float_t origin0[3]; // Origin of the generated parent particle (for GEANT tracking)
426 Float_t time0; // Time0 of the generated parent particle
427 Float_t pt, pl, ptot; // Transverse, logitudinal and total momenta of the parent particle
428 Float_t phi, theta; // Phi and theta spherical angles of the parent particle momentum
430 och[3]; // Momentum, polarisation and origin of the children particles from lujet
435 static TClonesArray *particles;
437 if(!particles) particles = new TClonesArray("TParticle",1000);
439 TDatabasePDG *pDataBase = TDatabasePDG::Instance();
443 // Calculating vertex position per event
444 for (j=0;j<3;j++) origin0[j]=fOrigin[j];
446 if(fVertexSmear==kPerEvent) {
448 for (j=0;j<3;j++) origin0[j]=fVertex[j];
454 // Generating fNpart particles
457 Int_t nGen = fNpart*ntimes;
462 Int_t iPart = fIpParaFunc(fRandom);
463 fChildWeight=(fDecayer->GetPartialBranchingRatio(iPart))*fParentWeight;
464 TParticlePDG *particle = pDataBase->GetParticle(iPart);
465 Float_t am = particle->Mass();
470 ty = TMath::TanH(fYPara->GetRandom());
473 if (fAnalog == kAnalog) {
474 pt=fPtPara->GetRandom();
478 pt=fPtMin+random[1]*(fPtMax-fPtMin);
480 wgtp=fParentWeight*fPtParaFunc(& ptd, &dummy);
481 wgtch=fChildWeight*fPtParaFunc(& ptd, &dummy);
483 xmt=sqrt(pt*pt+am*am);
484 if (TMath::Abs(ty)==1.) {
487 "Division by 0: Please check you rapidity range !");
492 //phi=fEvPlane; //align first particle of each event with event plane
494 double v2 = fV2Para->Eval(pt);
495 fdNdPhi->SetParameter(0,v2);
496 fdNdPhi->SetParameter(1,fEvPlane);
497 phi=fdNdPhi->GetRandom();
500 pl=xmt*ty/sqrt((1.-ty)*(1.+ty));
501 theta=TMath::ATan2(pt,pl);
503 if(theta<fThetaMin || theta>fThetaMax) continue;
504 ptot=TMath::Sqrt(pt*pt+pl*pl);
506 if(ptot<fPMin || ptot>fPMax) continue;
508 p[0]=pt*TMath::Cos(phi);
509 p[1]=pt*TMath::Sin(phi);
511 if(fVertexSmear==kPerTrack) {
515 fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
516 TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
519 time0 = fTimeOrigin + fOsigma[2]/TMath::Ccgs()*
520 TMath::Cos(2*random[0]*TMath::Pi())*
521 TMath::Sqrt(-2*TMath::Log(random[1]));
524 // Looking at fForceDecay :
525 // if fForceDecay != none Primary particle decays using
526 // AliPythia and children are tracked by GEANT
528 // if fForceDecay == none Primary particle is tracked by GEANT
529 // (In the latest, make sure that GEANT actually does all the decays you want)
531 Bool_t decayed = kFALSE;
534 if (fForceDecay != kNoDecay) {
535 // Using lujet to decay particle
536 Float_t energy=TMath::Sqrt(ptot*ptot+am*am);
537 TLorentzVector pmom(p[0], p[1], p[2], energy);
538 fDecayer->Decay(iPart,&pmom);
540 // select decay particles
541 Int_t np=fDecayer->ImportParticles(particles);
543 if(fForceConv) np=ForceGammaConversion(particles,np);
545 // for(int iPart=0; iPart<np; iPart++){
546 // TParticle *gamma = (TParticle *) particles->At(iPart);
547 // printf("%i %i:", iPart, gamma->GetPdgCode());
548 // printf("%i %i %i|",gamma->GetFirstMother(),gamma->GetFirstDaughter(),gamma->GetLastDaughter());
551 // Selecting GeometryAcceptance for particles fPdgCodeParticleforAcceptanceCut;
552 if (fGeometryAcceptance)
553 if (!CheckAcceptanceGeometry(np,particles)) continue;
555 Int_t* pFlag = new Int_t[np];
556 Int_t* pParent = new Int_t[np];
557 Int_t* pSelected = new Int_t[np];
558 Int_t* trackIt = new Int_t[np];
560 for (i=0; i<np; i++) {
568 TParticle* iparticle = 0;
570 for (i = 1; i<np ; i++) {
572 iparticle = (TParticle *) particles->At(i);
573 Int_t kf = iparticle->GetPdgCode();
574 Int_t ks = iparticle->GetStatusCode();
578 ipF = iparticle->GetFirstDaughter();
579 ipL = iparticle->GetLastDaughter();
580 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
584 // flag decay products of particles with long life-time (c tau > .3 mum)
587 // TParticlePDG *particle = pDataBase->GetParticle(kf);
589 Double_t lifeTime = fDecayer->GetLifetime(kf);
590 // Double_t mass = particle->Mass();
591 // Double_t width = particle->Width();
592 if (lifeTime > (Double_t) fMaxLifeTime) {
593 ipF = iparticle->GetFirstDaughter();
594 ipL = iparticle->GetLastDaughter();
595 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
604 if ((ChildSelected(TMath::Abs(kf)) || fForceDecay == kAll || fSelectAll) && trackIt[i])
607 pc[0]=iparticle->Px();
608 pc[1]=iparticle->Py();
609 pc[2]=iparticle->Pz();
610 Bool_t childok = KinematicSelection(iparticle, 1);
621 } // if child selection
623 } // decay particle loop
624 } // if decay products
627 if ((fCutOnChild && ncsel >0) || !fCutOnChild){
633 PushTrack(0, -1, iPart, p, origin0, polar, time0, kPPrimary, nt, wgtp, ((decayed)? 11 : 1));
641 for (i = 1; i < np; i++) {
643 TParticle* iparticle = (TParticle *) particles->At(i);
644 Int_t kf = iparticle->GetPdgCode();
645 Int_t ksc = iparticle->GetStatusCode();
646 Int_t jpa = iparticle->GetFirstMother()-1;
648 och[0] = origin0[0]+iparticle->Vx();
649 och[1] = origin0[1]+iparticle->Vy();
650 och[2] = origin0[2]+iparticle->Vz();
651 pc[0] = iparticle->Px();
652 pc[1] = iparticle->Py();
653 pc[2] = iparticle->Pz();
656 iparent = pParent[jpa];
661 PushTrack(fTrackIt * trackIt[i], iparent, kf,
663 time0 + iparticle->T(), kPDecay, nt, wgtch, ksc);
671 if (pFlag) delete[] pFlag;
672 if (pParent) delete[] pParent;
673 if (pSelected) delete[] pSelected;
674 if (trackIt) delete[] trackIt;
675 } // kinematic selection
676 else // nodecay option, so parent will be tracked by GEANT (pions, kaons, eta, omegas, baryons)
679 PushTrack(fTrackIt,-1,iPart,p,origin0,polar,time0,kPPrimary,nt,wgtp, 1);
687 SetHighWaterMark(nt);
689 AliGenEventHeader* header = new AliGenEventHeader("PARAM");
690 header->SetPrimaryVertex(fVertex);
691 header->SetInteractionTime(fTime);
692 header->SetNProduced(fNprimaries);
695 //____________________________________________________________________________________
696 Float_t AliGenParam::GetRelativeArea(Float_t ptMin, Float_t ptMax, Float_t yMin, Float_t yMax, Float_t phiMin, Float_t phiMax)
699 // Normalisation for selected kinematic region
701 #if ROOT_VERSION_CODE < ROOT_VERSION(5,99,0)
703 fPtPara->Integral(ptMin,ptMax,(Double_t *)0,1.e-6) / fPtPara->Integral( fPtPara->GetXmin(), fPtPara->GetXmax(),(Double_t *)0,1.e-6) *
704 fYPara->Integral(yMin,yMax,(Double_t *)0,1.e-6)/fYPara->Integral(fYPara->GetXmin(),fYPara->GetXmax(),(Double_t *)0,1.e-6) *
705 (phiMax-phiMin)/360.;
708 fPtPara->Integral(ptMin,ptMax,1.e-6) / fPtPara->Integral( fPtPara->GetXmin(), fPtPara->GetXmax(),1.e-6) *
709 fYPara->Integral(yMin,yMax,1.e-6)/fYPara->Integral(fYPara->GetXmin(),fYPara->GetXmax(),1.e-6) *
710 (phiMax-phiMin)/360.;
712 return TMath::Abs(ratio);
715 //____________________________________________________________________________________
717 void AliGenParam::Draw( const char * /*opt*/)
720 // Draw the pT and y Distributions
722 TCanvas *c0 = new TCanvas("c0","Canvas 0",400,10,600,700);
726 fPtPara->GetHistogram()->SetXTitle("p_{T} (GeV)");
729 fYPara->GetHistogram()->SetXTitle("y");