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 <TDatabasePDG.h>
28 #include <TLorentzVector.h>
30 #include <TParticle.h>
31 #include <TParticlePDG.h>
33 #include <TVirtualMC.h>
35 #include "AliDecayer.h"
36 #include "AliGenMUONlib.h"
37 #include "AliGenParam.h"
43 //------------------------------------------------------------
47 <img src="picts/AliGenParam.gif">
51 //____________________________________________________________
52 AliGenParam::AliGenParam():
67 // Default constructor
69 //____________________________________________________________
70 AliGenParam::AliGenParam(Int_t npart, AliGenLib * Library, Int_t param, char* tname)
72 fPtParaFunc(Library->GetPt(param, tname)),
73 fYParaFunc (Library->GetY (param, tname)),
74 fIpParaFunc(Library->GetIp(param, tname)),
86 // Constructor using number of particles parameterisation id and library
88 fTitle= "Particle Generator using pT and y parameterisation";
92 //____________________________________________________________
93 AliGenParam::AliGenParam(Int_t npart, Int_t param, const char* tname, const char* name):
109 // Constructor using parameterisation id and number of particles
112 fTitle= "Particle Generator using pT and y parameterisation";
114 AliGenLib* pLibrary = new AliGenMUONlib();
115 fPtParaFunc = pLibrary->GetPt(param, tname);
116 fYParaFunc = pLibrary->GetY (param, tname);
117 fIpParaFunc = pLibrary->GetIp(param, tname);
121 for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
124 SetChildMomentumRange();
127 SetChildThetaRange();
129 //____________________________________________________________
131 AliGenParam::AliGenParam(Int_t npart, Int_t param,
132 Double_t (*PtPara) (Double_t*, Double_t*),
133 Double_t (*YPara ) (Double_t* ,Double_t*),
134 Int_t (*IpPara) (TRandom *))
152 // Gines Martinez 1/10/99
154 fTitle= "Particle Generator using pT and y parameterisation";
158 for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
161 SetChildMomentumRange();
164 SetChildThetaRange();
168 AliGenParam::AliGenParam(const AliGenParam & Param)
188 //____________________________________________________________
189 AliGenParam::~AliGenParam()
196 //____________________________________________________________
197 void AliGenParam::Init()
201 if (gMC) fDecayer = gMC->GetDecayer();
204 <img src="picts/AliGenParam.gif">
208 sprintf(name, "pt-parameterisation for %s", GetName());
210 if (fPtPara) fPtPara->Delete();
211 fPtPara = new TF1(name, fPtParaFunc, fPtMin, fPtMax,0);
212 gROOT->GetListOfFunctions()->Remove(fPtPara);
213 // Set representation precision to 10 MeV
214 Int_t npx= Int_t((fPtMax - fPtMin) / fDeltaPt);
216 fPtPara->SetNpx(npx);
218 sprintf(name, "y-parameterisation for %s", GetName());
219 if (fYPara) fYPara->Delete();
220 fYPara = new TF1(name, fYParaFunc, fYMin, fYMax, 0);
221 gROOT->GetListOfFunctions()->Remove(fYPara);
224 sprintf(name, "pt-for-%s", GetName());
225 TF1 ptPara(name ,fPtParaFunc, 0, 15, 0);
226 sprintf(name, "y-for-%s", GetName());
227 TF1 yPara(name, fYParaFunc, -6, 6, 0);
234 fdNdy0=fYParaFunc(&y1,&y2);
236 // Integral over generation region
237 Float_t intYS = yPara.Integral(fYMin, fYMax,(Double_t*) 0x0,1.e-6);
238 Float_t intPt0 = ptPara.Integral(0,15,(Double_t *) 0x0,1.e-6);
239 Float_t intPtS = ptPara.Integral(fPtMin,fPtMax,(Double_t*) 0x0,1.e-6);
240 Float_t phiWgt=(fPhiMax-fPhiMin)/2./TMath::Pi();
241 if (fAnalog == kAnalog) {
242 fYWgt = intYS/fdNdy0;
243 fPtWgt = intPtS/intPt0;
244 fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
246 fYWgt = intYS/fdNdy0;
247 fPtWgt = (fPtMax-fPtMin)/intPt0;
248 fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
251 // particle decay related initialization
252 fDecayer->SetForceDecay(fForceDecay);
259 //____________________________________________________________
260 void AliGenParam::Generate()
263 // Generate 'npart' of light and heavy mesons (J/Psi, upsilon or phi, Pion,
264 // Kaons, Etas, Omegas) and Baryons (proton, antiprotons, neutrons and
265 // antineutrons in the the desired theta, phi and momentum windows;
266 // Gaussian smearing on the vertex is done if selected.
267 // The decay of heavy mesons is done using lujet,
268 // and the childern particle are tracked by GEANT
269 // However, light mesons are directly tracked by GEANT
270 // setting fForceDecay = nodecay (SetForceDecay(nodecay))
273 // Reinitialize decayer
274 fDecayer->SetForceDecay(fForceDecay);
278 Float_t polar[3]= {0,0,0}; // Polarisation of the parent particle (for GEANT tracking)
279 Float_t origin0[3]; // Origin of the generated parent particle (for GEANT tracking)
280 Float_t pt, pl, ptot; // Transverse, logitudinal and total momenta of the parent particle
281 Float_t phi, theta; // Phi and theta spherical angles of the parent particle momentum
283 och[3]; // Momentum, polarisation and origin of the children particles from lujet
288 static TClonesArray *particles;
290 if(!particles) particles = new TClonesArray("TParticle",1000);
292 TDatabasePDG *pDataBase = TDatabasePDG::Instance();
296 // Calculating vertex position per event
297 for (j=0;j<3;j++) origin0[j]=fOrigin[j];
298 if(fVertexSmear==kPerEvent) {
300 for (j=0;j<3;j++) origin0[j]=fVertex[j];
305 // Generating fNpart particles
310 Int_t iPart = fIpParaFunc(fRandom);
311 fChildWeight=(fDecayer->GetPartialBranchingRatio(iPart))*fParentWeight;
312 TParticlePDG *particle = pDataBase->GetParticle(iPart);
313 Float_t am = particle->Mass();
318 phi=fPhiMin+random[0]*(fPhiMax-fPhiMin);
321 ty = TMath::TanH(fYPara->GetRandom());
324 if (fAnalog == kAnalog) {
325 pt=fPtPara->GetRandom();
329 pt=fPtMin+random[1]*(fPtMax-fPtMin);
331 wgtp=fParentWeight*fPtParaFunc(& ptd, &dummy);
332 wgtch=fChildWeight*fPtParaFunc(& ptd, &dummy);
334 xmt=sqrt(pt*pt+am*am);
335 if (TMath::Abs(ty)==1.) {
338 "Division by 0: Please check you rapidity range !");
341 pl=xmt*ty/sqrt(1.-ty*ty);
342 theta=TMath::ATan2(pt,pl);
344 if(theta<fThetaMin || theta>fThetaMax) continue;
345 ptot=TMath::Sqrt(pt*pt+pl*pl);
347 if(ptot<fPMin || ptot>fPMax) continue;
349 p[0]=pt*TMath::Cos(phi);
350 p[1]=pt*TMath::Sin(phi);
352 if(fVertexSmear==kPerTrack) {
356 fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
357 TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
361 // Looking at fForceDecay :
362 // if fForceDecay != none Primary particle decays using
363 // AliPythia and children are tracked by GEANT
365 // if fForceDecay == none Primary particle is tracked by GEANT
366 // (In the latest, make sure that GEANT actually does all the decays you want)
369 if (fForceDecay != kNoDecay) {
370 // Using lujet to decay particle
371 Float_t energy=TMath::Sqrt(ptot*ptot+am*am);
372 TLorentzVector pmom(p[0], p[1], p[2], energy);
373 fDecayer->Decay(iPart,&pmom);
375 // select decay particles
376 Int_t np=fDecayer->ImportParticles(particles);
378 // Selecting GeometryAcceptance for particles fPdgCodeParticleforAcceptanceCut;
379 if (fGeometryAcceptance)
380 if (!CheckAcceptanceGeometry(np,particles)) continue;
382 Int_t* pFlag = new Int_t[np];
383 Int_t* pParent = new Int_t[np];
384 Int_t* pSelected = new Int_t[np];
385 Int_t* trackIt = new Int_t[np];
387 for (i=0; i<np; i++) {
394 TParticle* iparticle = (TParticle *) particles->At(0);
396 for (i = 1; i<np ; i++) {
398 iparticle = (TParticle *) particles->At(i);
399 Int_t kf = iparticle->GetPdgCode();
400 Int_t ks = iparticle->GetStatusCode();
404 ipF = iparticle->GetFirstDaughter();
405 ipL = iparticle->GetLastDaughter();
406 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
410 // flag decay products of particles with long life-time (c tau > .3 mum)
413 // TParticlePDG *particle = pDataBase->GetParticle(kf);
415 Double_t lifeTime = fDecayer->GetLifetime(kf);
416 // Double_t mass = particle->Mass();
417 // Double_t width = particle->Width();
418 if (lifeTime > (Double_t) fMaxLifeTime) {
419 ipF = iparticle->GetFirstDaughter();
420 ipL = iparticle->GetLastDaughter();
421 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
430 if (ChildSelected(TMath::Abs(kf)) || fForceDecay == kAll && trackIt[i])
433 pc[0]=iparticle->Px();
434 pc[1]=iparticle->Py();
435 pc[2]=iparticle->Pz();
436 Bool_t childok = KinematicSelection(iparticle, 1);
447 } // if child selection
449 } // decay particle loop
450 } // if decay products
453 if ((fCutOnChild && ncsel >0) || !fCutOnChild){
457 PushTrack(0, -1, iPart, p, origin0, polar, 0, kPPrimary, nt, wgtp);
463 for (i = 1; i < np; i++) {
465 TParticle* iparticle = (TParticle *) particles->At(i);
466 Int_t kf = iparticle->GetPdgCode();
467 Int_t ipa = iparticle->GetFirstMother()-1;
469 och[0] = origin0[0]+iparticle->Vx()/10;
470 och[1] = origin0[1]+iparticle->Vy()/10;
471 och[2] = origin0[2]+iparticle->Vz()/10;
472 pc[0] = iparticle->Px();
473 pc[1] = iparticle->Py();
474 pc[2] = iparticle->Pz();
477 iparent = pParent[ipa];
482 PushTrack(fTrackIt*trackIt[i], iparent, kf,
484 0, kPDecay, nt, wgtch);
491 if (pFlag) delete[] pFlag;
492 if (pParent) delete[] pParent;
493 if (pSelected) delete[] pSelected;
494 if (trackIt) delete[] trackIt;
495 } // kinematic selection
496 else // nodecay option, so parent will be tracked by GEANT (pions, kaons, eta, omegas, baryons)
499 PushTrack(fTrackIt,-1,iPart,p,origin0,polar,0,kPPrimary,nt,wgtp);
505 SetHighWaterMark(nt);
507 //____________________________________________________________________________________
508 Float_t AliGenParam::GetRelativeArea(Float_t ptMin, Float_t ptMax, Float_t yMin, Float_t yMax, Float_t phiMin, Float_t phiMax)
511 // Normalisation for selected kinematic region
514 fPtPara->Integral(ptMin,ptMax,(Double_t *)0,1.e-6) / fPtPara->Integral( fPtPara->GetXmin(), fPtPara->GetXmax(),(Double_t *)0,1.e-6) *
515 fYPara->Integral(yMin,yMax,(Double_t *)0,1.e-6)/fYPara->Integral(fYPara->GetXmin(),fYPara->GetXmax(),(Double_t *)0,1.e-6) *
516 (phiMax-phiMin)/360.;
517 return TMath::Abs(ratio);
520 //____________________________________________________________________________________
522 void AliGenParam::Draw( const char * /*opt*/)
525 // Draw the pT and y Distributions
527 TCanvas *c0 = new TCanvas("c0","Canvas 0",400,10,600,700);
531 fPtPara->GetHistogram()->SetXTitle("p_{T} (GeV)");
534 fYPara->GetHistogram()->SetXTitle("y");
537 AliGenParam& AliGenParam::operator=(const AliGenParam& rhs)
539 // Assignment operator