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
24 #include "AliGenParam.h"
25 #include "AliDecayer.h"
26 #include "AliGenMUONlib.h"
28 #include <TParticle.h>
29 #include <TParticlePDG.h>
30 #include <TDatabasePDG.h>
31 #include <TLorentzVector.h>
32 #include <TVirtualMC.h>
42 //------------------------------------------------------------
46 <img src="picts/AliGenParam.gif">
50 //____________________________________________________________
51 AliGenParam::AliGenParam()
53 // Deafault constructor
63 //____________________________________________________________
64 AliGenParam::AliGenParam(Int_t npart, AliGenLib * Library, Int_t param, char* tname):AliGenMC(npart)
66 // Constructor using number of particles parameterisation id and library
68 fTitle= "Particle Generator using pT and y parameterisation";
70 fPtParaFunc = Library->GetPt(param, tname);
71 fYParaFunc = Library->GetY (param, tname);
72 fIpParaFunc = Library->GetIp(param, tname);
81 //____________________________________________________________
82 AliGenParam::AliGenParam(Int_t npart, Int_t param, const char* tname, const char* name):AliGenMC(npart)
84 // Constructor using parameterisation id and number of particles
87 fTitle= "Particle Generator using pT and y parameterisation";
89 AliGenLib* pLibrary = new AliGenMUONlib();
91 fPtParaFunc = pLibrary->GetPt(param, tname);
92 fYParaFunc = pLibrary->GetY (param, tname);
93 fIpParaFunc = pLibrary->GetIp(param, tname);
100 for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
103 SetChildMomentumRange();
106 SetChildThetaRange();
109 //____________________________________________________________
111 AliGenParam::AliGenParam(Int_t npart, Int_t param,
112 Double_t (*PtPara) (Double_t*, Double_t*),
113 Double_t (*YPara ) (Double_t* ,Double_t*),
114 Int_t (*IpPara) (TRandom *))
118 // Gines Martinez 1/10/99
120 fTitle= "Particle Generator using pT and y parameterisation";
122 fPtParaFunc = PtPara;
124 fIpParaFunc = IpPara;
131 for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
134 SetChildMomentumRange();
137 SetChildThetaRange();
142 AliGenParam::AliGenParam(const AliGenParam & Param)
149 //____________________________________________________________
150 AliGenParam::~AliGenParam()
157 //____________________________________________________________
158 void AliGenParam::Init()
162 if (gMC) fDecayer = gMC->GetDecayer();
165 <img src="picts/AliGenParam.gif">
169 sprintf(name, "pt-parameterisation for %s", GetName());
171 if (fPtPara) fPtPara->Delete();
172 fPtPara = new TF1(name, fPtParaFunc, fPtMin, fPtMax,0);
173 // Set representation precision to 10 MeV
174 Int_t npx= Int_t((fPtMax - fPtMin) / fDeltaPt);
176 fPtPara->SetNpx(npx);
178 sprintf(name, "y-parameterisation for %s", GetName());
179 if (fYPara) fYPara->Delete();
180 fYPara = new TF1(name, fYParaFunc, fYMin, fYMax, 0);
182 sprintf(name, "pt-for-%s", GetName());
183 TF1 ptPara(name ,fPtParaFunc, 0, 15, 0);
184 sprintf(name, "y-for-%s", GetName());
185 TF1 yPara(name, fYParaFunc, -6, 6, 0);
192 fdNdy0=fYParaFunc(&y1,&y2);
194 // Integral over generation region
195 Float_t intYS = yPara.Integral(fYMin, fYMax);
196 Float_t intPt0 = ptPara.Integral(0,15);
197 Float_t intPtS = ptPara.Integral(fPtMin,fPtMax);
198 Float_t phiWgt=(fPhiMax-fPhiMin)/2./TMath::Pi();
199 if (fAnalog == kAnalog) {
200 fYWgt = intYS/fdNdy0;
201 fPtWgt = intPtS/intPt0;
202 fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
204 fYWgt = intYS/fdNdy0;
205 fPtWgt = (fPtMax-fPtMin)/intPt0;
206 fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
209 // particle decay related initialization
210 fDecayer->SetForceDecay(fForceDecay);
217 //____________________________________________________________
218 void AliGenParam::Generate()
221 // Generate 'npart' of light and heavy mesons (J/Psi, upsilon or phi, Pion,
222 // Kaons, Etas, Omegas) and Baryons (proton, antiprotons, neutrons and
223 // antineutrons in the the desired theta, phi and momentum windows;
224 // Gaussian smearing on the vertex is done if selected.
225 // The decay of heavy mesons is done using lujet,
226 // and the childern particle are tracked by GEANT
227 // However, light mesons are directly tracked by GEANT
228 // setting fForceDecay = nodecay (SetForceDecay(nodecay))
231 // Reinitialize decayer
232 fDecayer->SetForceDecay(fForceDecay);
236 Float_t polar[3]= {0,0,0}; // Polarisation of the parent particle (for GEANT tracking)
237 Float_t origin0[3]; // Origin of the generated parent particle (for GEANT tracking)
238 Float_t pt, pl, ptot; // Transverse, logitudinal and total momenta of the parent particle
239 Float_t phi, theta; // Phi and theta spherical angles of the parent particle momentum
241 och[3]; // Momentum, polarisation and origin of the children particles from lujet
246 static TClonesArray *particles;
248 if(!particles) particles = new TClonesArray("TParticle",1000);
250 TDatabasePDG *pDataBase = TDatabasePDG::Instance();
254 // Calculating vertex position per event
255 for (j=0;j<3;j++) origin0[j]=fOrigin[j];
256 if(fVertexSmear==kPerEvent) {
258 for (j=0;j<3;j++) origin0[j]=fVertex[j];
263 // Generating fNpart particles
268 Int_t iPart = fIpParaFunc(fRandom);
269 fChildWeight=(fDecayer->GetPartialBranchingRatio(iPart))*fParentWeight;
270 TParticlePDG *particle = pDataBase->GetParticle(iPart);
271 Float_t am = particle->Mass();
276 phi=fPhiMin+random[0]*(fPhiMax-fPhiMin);
279 ty = TMath::TanH(fYPara->GetRandom());
282 if (fAnalog == kAnalog) {
283 pt=fPtPara->GetRandom();
287 pt=fPtMin+random[1]*(fPtMax-fPtMin);
289 wgtp=fParentWeight*fPtParaFunc(& ptd, &dummy);
290 wgtch=fChildWeight*fPtParaFunc(& ptd, &dummy);
292 xmt=sqrt(pt*pt+am*am);
293 if (TMath::Abs(ty)==1.) {
296 "Division by 0: Please check you rapidity range !");
299 pl=xmt*ty/sqrt(1.-ty*ty);
300 theta=TMath::ATan2(pt,pl);
302 if(theta<fThetaMin || theta>fThetaMax) continue;
303 ptot=TMath::Sqrt(pt*pt+pl*pl);
305 if(ptot<fPMin || ptot>fPMax) continue;
307 p[0]=pt*TMath::Cos(phi);
308 p[1]=pt*TMath::Sin(phi);
310 if(fVertexSmear==kPerTrack) {
314 fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
315 TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
319 // Looking at fForceDecay :
320 // if fForceDecay != none Primary particle decays using
321 // AliPythia and children are tracked by GEANT
323 // if fForceDecay == none Primary particle is tracked by GEANT
324 // (In the latest, make sure that GEANT actually does all the decays you want)
327 if (fForceDecay != kNoDecay) {
328 // Using lujet to decay particle
329 Float_t energy=TMath::Sqrt(ptot*ptot+am*am);
330 TLorentzVector pmom(p[0], p[1], p[2], energy);
331 fDecayer->Decay(iPart,&pmom);
333 // select decay particles
334 Int_t np=fDecayer->ImportParticles(particles);
336 // Selecting GeometryAcceptance for particles fPdgCodeParticleforAcceptanceCut;
337 if (fGeometryAcceptance)
338 if (!CheckAcceptanceGeometry(np,particles)) continue;
340 Int_t* pFlag = new Int_t[np];
341 Int_t* pParent = new Int_t[np];
342 Int_t* pSelected = new Int_t[np];
343 Int_t* trackIt = new Int_t[np];
345 for (i=0; i<np; i++) {
352 TParticle* iparticle = (TParticle *) particles->At(0);
354 for (i = 1; i<np ; i++) {
356 iparticle = (TParticle *) particles->At(i);
357 Int_t kf = iparticle->GetPdgCode();
358 Int_t ks = iparticle->GetStatusCode();
362 ipF = iparticle->GetFirstDaughter();
363 ipL = iparticle->GetLastDaughter();
364 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
368 // flag decay products of particles with long life-time (c tau > .3 mum)
371 // TParticlePDG *particle = pDataBase->GetParticle(kf);
373 Double_t lifeTime = fDecayer->GetLifetime(kf);
374 // Double_t mass = particle->Mass();
375 // Double_t width = particle->Width();
376 if (lifeTime > (Double_t) fMaxLifeTime) {
377 ipF = iparticle->GetFirstDaughter();
378 ipL = iparticle->GetLastDaughter();
379 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
388 if (ChildSelected(TMath::Abs(kf)) || fForceDecay == kAll && trackIt[i])
391 pc[0]=iparticle->Px();
392 pc[1]=iparticle->Py();
393 pc[2]=iparticle->Pz();
394 Bool_t childok = KinematicSelection(iparticle, 1);
405 } // if child selection
407 } // decay particle loop
408 } // if decay products
411 if ((fCutOnChild && ncsel >0) || !fCutOnChild){
415 PushTrack(0, -1, iPart, p, origin0, polar, 0, kPPrimary, nt, wgtp);
421 for (i = 1; i < np; i++) {
423 TParticle* iparticle = (TParticle *) particles->At(i);
424 Int_t kf = iparticle->GetPdgCode();
425 Int_t ipa = iparticle->GetFirstMother()-1;
427 och[0] = origin0[0]+iparticle->Vx()/10;
428 och[1] = origin0[1]+iparticle->Vy()/10;
429 och[2] = origin0[2]+iparticle->Vz()/10;
430 pc[0] = iparticle->Px();
431 pc[1] = iparticle->Py();
432 pc[2] = iparticle->Pz();
435 iparent = pParent[ipa];
440 PushTrack(fTrackIt*trackIt[i], iparent, kf,
442 0, kPDecay, nt, wgtch);
449 if (pFlag) delete[] pFlag;
450 if (pParent) delete[] pParent;
451 if (pSelected) delete[] pSelected;
452 if (trackIt) delete[] trackIt;
453 } // kinematic selection
454 else // nodecay option, so parent will be tracked by GEANT (pions, kaons, eta, omegas, baryons)
457 PushTrack(fTrackIt,-1,iPart,p,origin0,polar,0,kPPrimary,nt,wgtp);
463 SetHighWaterMark(nt);
465 //____________________________________________________________________________________
466 Float_t AliGenParam::GetRelativeArea(Float_t ptMin, Float_t ptMax, Float_t yMin, Float_t yMax, Float_t phiMin, Float_t phiMax)
469 // Normalisation for selected kinematic region
472 fPtPara->Integral(ptMin,ptMax) / fPtPara->Integral( fPtPara->GetXmin(), fPtPara->GetXmax()) *
473 fYPara->Integral(yMin,yMax)/fYPara->Integral(fYPara->GetXmin(),fYPara->GetXmax()) *
474 (phiMax-phiMin)/360.;
475 return TMath::Abs(ratio);
478 //____________________________________________________________________________________
480 void AliGenParam::Draw( const char * /*opt*/)
483 // Draw the pT and y Distributions
485 TCanvas *c0 = new TCanvas("c0","Canvas 0",400,10,600,700);
489 fPtPara->GetHistogram()->SetXTitle("p_{T} (GeV)");
492 fYPara->GetHistogram()->SetXTitle("y");
495 AliGenParam& AliGenParam::operator=(const AliGenParam& rhs)
497 // Assignment operator