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()
54 // Deafault constructor
64 //____________________________________________________________
65 AliGenParam::AliGenParam(Int_t npart, AliGenLib * Library, Int_t param, char* tname):AliGenMC(npart)
67 // Constructor using number of particles parameterisation id and library
69 fTitle= "Particle Generator using pT and y parameterisation";
71 fPtParaFunc = Library->GetPt(param, tname);
72 fYParaFunc = Library->GetY (param, tname);
73 fIpParaFunc = Library->GetIp(param, tname);
82 //____________________________________________________________
83 AliGenParam::AliGenParam(Int_t npart, Int_t param, const char* tname, const char* name):AliGenMC(npart)
85 // Constructor using parameterisation id and number of particles
88 fTitle= "Particle Generator using pT and y parameterisation";
90 AliGenLib* pLibrary = new AliGenMUONlib();
92 fPtParaFunc = pLibrary->GetPt(param, tname);
93 fYParaFunc = pLibrary->GetY (param, tname);
94 fIpParaFunc = pLibrary->GetIp(param, tname);
101 for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
104 SetChildMomentumRange();
107 SetChildThetaRange();
110 //____________________________________________________________
112 AliGenParam::AliGenParam(Int_t npart, Int_t param,
113 Double_t (*PtPara) (Double_t*, Double_t*),
114 Double_t (*YPara ) (Double_t* ,Double_t*),
115 Int_t (*IpPara) (TRandom *))
119 // Gines Martinez 1/10/99
121 fTitle= "Particle Generator using pT and y parameterisation";
123 fPtParaFunc = PtPara;
125 fIpParaFunc = IpPara;
132 for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
135 SetChildMomentumRange();
138 SetChildThetaRange();
143 AliGenParam::AliGenParam(const AliGenParam & Param)
150 //____________________________________________________________
151 AliGenParam::~AliGenParam()
158 //____________________________________________________________
159 void AliGenParam::Init()
163 if (gMC) fDecayer = gMC->GetDecayer();
166 <img src="picts/AliGenParam.gif">
170 sprintf(name, "pt-parameterisation for %s", GetName());
172 if (fPtPara) fPtPara->Delete();
173 fPtPara = new TF1(name, fPtParaFunc, fPtMin, fPtMax,0);
174 gROOT->GetListOfFunctions()->Remove(fPtPara);
175 // Set representation precision to 10 MeV
176 Int_t npx= Int_t((fPtMax - fPtMin) / fDeltaPt);
178 fPtPara->SetNpx(npx);
180 sprintf(name, "y-parameterisation for %s", GetName());
181 if (fYPara) fYPara->Delete();
182 fYPara = new TF1(name, fYParaFunc, fYMin, fYMax, 0);
183 gROOT->GetListOfFunctions()->Remove(fYPara);
186 sprintf(name, "pt-for-%s", GetName());
187 TF1 ptPara(name ,fPtParaFunc, 0, 15, 0);
188 sprintf(name, "y-for-%s", GetName());
189 TF1 yPara(name, fYParaFunc, -6, 6, 0);
196 fdNdy0=fYParaFunc(&y1,&y2);
198 // Integral over generation region
199 Float_t intYS = yPara.Integral(fYMin, fYMax,(Double_t*) 0x0,1.e-6);
200 Float_t intPt0 = ptPara.Integral(0,15,(Double_t *) 0x0,1.e-6);
201 Float_t intPtS = ptPara.Integral(fPtMin,fPtMax,(Double_t*) 0x0,1.e-6);
202 Float_t phiWgt=(fPhiMax-fPhiMin)/2./TMath::Pi();
203 if (fAnalog == kAnalog) {
204 fYWgt = intYS/fdNdy0;
205 fPtWgt = intPtS/intPt0;
206 fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
208 fYWgt = intYS/fdNdy0;
209 fPtWgt = (fPtMax-fPtMin)/intPt0;
210 fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
213 // particle decay related initialization
214 fDecayer->SetForceDecay(fForceDecay);
221 //____________________________________________________________
222 void AliGenParam::Generate()
225 // Generate 'npart' of light and heavy mesons (J/Psi, upsilon or phi, Pion,
226 // Kaons, Etas, Omegas) and Baryons (proton, antiprotons, neutrons and
227 // antineutrons in the the desired theta, phi and momentum windows;
228 // Gaussian smearing on the vertex is done if selected.
229 // The decay of heavy mesons is done using lujet,
230 // and the childern particle are tracked by GEANT
231 // However, light mesons are directly tracked by GEANT
232 // setting fForceDecay = nodecay (SetForceDecay(nodecay))
235 // Reinitialize decayer
236 fDecayer->SetForceDecay(fForceDecay);
240 Float_t polar[3]= {0,0,0}; // Polarisation of the parent particle (for GEANT tracking)
241 Float_t origin0[3]; // Origin of the generated parent particle (for GEANT tracking)
242 Float_t pt, pl, ptot; // Transverse, logitudinal and total momenta of the parent particle
243 Float_t phi, theta; // Phi and theta spherical angles of the parent particle momentum
245 och[3]; // Momentum, polarisation and origin of the children particles from lujet
250 static TClonesArray *particles;
252 if(!particles) particles = new TClonesArray("TParticle",1000);
254 TDatabasePDG *pDataBase = TDatabasePDG::Instance();
258 // Calculating vertex position per event
259 for (j=0;j<3;j++) origin0[j]=fOrigin[j];
260 if(fVertexSmear==kPerEvent) {
262 for (j=0;j<3;j++) origin0[j]=fVertex[j];
267 // Generating fNpart particles
272 Int_t iPart = fIpParaFunc(fRandom);
273 fChildWeight=(fDecayer->GetPartialBranchingRatio(iPart))*fParentWeight;
274 TParticlePDG *particle = pDataBase->GetParticle(iPart);
275 Float_t am = particle->Mass();
280 phi=fPhiMin+random[0]*(fPhiMax-fPhiMin);
283 ty = TMath::TanH(fYPara->GetRandom());
286 if (fAnalog == kAnalog) {
287 pt=fPtPara->GetRandom();
291 pt=fPtMin+random[1]*(fPtMax-fPtMin);
293 wgtp=fParentWeight*fPtParaFunc(& ptd, &dummy);
294 wgtch=fChildWeight*fPtParaFunc(& ptd, &dummy);
296 xmt=sqrt(pt*pt+am*am);
297 if (TMath::Abs(ty)==1.) {
300 "Division by 0: Please check you rapidity range !");
303 pl=xmt*ty/sqrt(1.-ty*ty);
304 theta=TMath::ATan2(pt,pl);
306 if(theta<fThetaMin || theta>fThetaMax) continue;
307 ptot=TMath::Sqrt(pt*pt+pl*pl);
309 if(ptot<fPMin || ptot>fPMax) continue;
311 p[0]=pt*TMath::Cos(phi);
312 p[1]=pt*TMath::Sin(phi);
314 if(fVertexSmear==kPerTrack) {
318 fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
319 TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
323 // Looking at fForceDecay :
324 // if fForceDecay != none Primary particle decays using
325 // AliPythia and children are tracked by GEANT
327 // if fForceDecay == none Primary particle is tracked by GEANT
328 // (In the latest, make sure that GEANT actually does all the decays you want)
331 if (fForceDecay != kNoDecay) {
332 // Using lujet to decay particle
333 Float_t energy=TMath::Sqrt(ptot*ptot+am*am);
334 TLorentzVector pmom(p[0], p[1], p[2], energy);
335 fDecayer->Decay(iPart,&pmom);
337 // select decay particles
338 Int_t np=fDecayer->ImportParticles(particles);
340 // Selecting GeometryAcceptance for particles fPdgCodeParticleforAcceptanceCut;
341 if (fGeometryAcceptance)
342 if (!CheckAcceptanceGeometry(np,particles)) continue;
344 Int_t* pFlag = new Int_t[np];
345 Int_t* pParent = new Int_t[np];
346 Int_t* pSelected = new Int_t[np];
347 Int_t* trackIt = new Int_t[np];
349 for (i=0; i<np; i++) {
356 TParticle* iparticle = (TParticle *) particles->At(0);
358 for (i = 1; i<np ; i++) {
360 iparticle = (TParticle *) particles->At(i);
361 Int_t kf = iparticle->GetPdgCode();
362 Int_t ks = iparticle->GetStatusCode();
366 ipF = iparticle->GetFirstDaughter();
367 ipL = iparticle->GetLastDaughter();
368 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
372 // flag decay products of particles with long life-time (c tau > .3 mum)
375 // TParticlePDG *particle = pDataBase->GetParticle(kf);
377 Double_t lifeTime = fDecayer->GetLifetime(kf);
378 // Double_t mass = particle->Mass();
379 // Double_t width = particle->Width();
380 if (lifeTime > (Double_t) fMaxLifeTime) {
381 ipF = iparticle->GetFirstDaughter();
382 ipL = iparticle->GetLastDaughter();
383 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
392 if (ChildSelected(TMath::Abs(kf)) || fForceDecay == kAll && trackIt[i])
395 pc[0]=iparticle->Px();
396 pc[1]=iparticle->Py();
397 pc[2]=iparticle->Pz();
398 Bool_t childok = KinematicSelection(iparticle, 1);
409 } // if child selection
411 } // decay particle loop
412 } // if decay products
415 if ((fCutOnChild && ncsel >0) || !fCutOnChild){
419 PushTrack(0, -1, iPart, p, origin0, polar, 0, kPPrimary, nt, wgtp);
425 for (i = 1; i < np; i++) {
427 TParticle* iparticle = (TParticle *) particles->At(i);
428 Int_t kf = iparticle->GetPdgCode();
429 Int_t ipa = iparticle->GetFirstMother()-1;
431 och[0] = origin0[0]+iparticle->Vx()/10;
432 och[1] = origin0[1]+iparticle->Vy()/10;
433 och[2] = origin0[2]+iparticle->Vz()/10;
434 pc[0] = iparticle->Px();
435 pc[1] = iparticle->Py();
436 pc[2] = iparticle->Pz();
439 iparent = pParent[ipa];
444 PushTrack(fTrackIt*trackIt[i], iparent, kf,
446 0, kPDecay, nt, wgtch);
453 if (pFlag) delete[] pFlag;
454 if (pParent) delete[] pParent;
455 if (pSelected) delete[] pSelected;
456 if (trackIt) delete[] trackIt;
457 } // kinematic selection
458 else // nodecay option, so parent will be tracked by GEANT (pions, kaons, eta, omegas, baryons)
461 PushTrack(fTrackIt,-1,iPart,p,origin0,polar,0,kPPrimary,nt,wgtp);
467 SetHighWaterMark(nt);
469 //____________________________________________________________________________________
470 Float_t AliGenParam::GetRelativeArea(Float_t ptMin, Float_t ptMax, Float_t yMin, Float_t yMax, Float_t phiMin, Float_t phiMax)
473 // Normalisation for selected kinematic region
476 fPtPara->Integral(ptMin,ptMax,(Double_t *)0,1.e-6) / fPtPara->Integral( fPtPara->GetXmin(), fPtPara->GetXmax(),(Double_t *)0,1.e-6) *
477 fYPara->Integral(yMin,yMax,(Double_t *)0,1.e-6)/fYPara->Integral(fYPara->GetXmin(),fYPara->GetXmax(),(Double_t *)0,1.e-6) *
478 (phiMax-phiMin)/360.;
479 return TMath::Abs(ratio);
482 //____________________________________________________________________________________
484 void AliGenParam::Draw( const char * /*opt*/)
487 // Draw the pT and y Distributions
489 TCanvas *c0 = new TCanvas("c0","Canvas 0",400,10,600,700);
493 fPtPara->GetHistogram()->SetXTitle("p_{T} (GeV)");
496 fYPara->GetHistogram()->SetXTitle("y");
499 AliGenParam& AliGenParam::operator=(const AliGenParam& rhs)
501 // Assignment operator