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();
167 //____________________________________________________________
168 AliGenParam::~AliGenParam()
175 //____________________________________________________________
176 void AliGenParam::Init()
180 if (gMC) fDecayer = gMC->GetDecayer();
183 <img src="picts/AliGenParam.gif">
187 sprintf(name, "pt-parameterisation for %s", GetName());
189 if (fPtPara) fPtPara->Delete();
190 fPtPara = new TF1(name, fPtParaFunc, fPtMin, fPtMax,0);
191 gROOT->GetListOfFunctions()->Remove(fPtPara);
192 // Set representation precision to 10 MeV
193 Int_t npx= Int_t((fPtMax - fPtMin) / fDeltaPt);
195 fPtPara->SetNpx(npx);
197 sprintf(name, "y-parameterisation for %s", GetName());
198 if (fYPara) fYPara->Delete();
199 fYPara = new TF1(name, fYParaFunc, fYMin, fYMax, 0);
200 gROOT->GetListOfFunctions()->Remove(fYPara);
203 sprintf(name, "pt-for-%s", GetName());
204 TF1 ptPara(name ,fPtParaFunc, 0, 15, 0);
205 sprintf(name, "y-for-%s", GetName());
206 TF1 yPara(name, fYParaFunc, -6, 6, 0);
213 fdNdy0=fYParaFunc(&y1,&y2);
215 // Integral over generation region
216 Float_t intYS = yPara.Integral(fYMin, fYMax,(Double_t*) 0x0,1.e-6);
217 Float_t intPt0 = ptPara.Integral(0,15,(Double_t *) 0x0,1.e-6);
218 Float_t intPtS = ptPara.Integral(fPtMin,fPtMax,(Double_t*) 0x0,1.e-6);
219 Float_t phiWgt=(fPhiMax-fPhiMin)/2./TMath::Pi();
220 if (fAnalog == kAnalog) {
221 fYWgt = intYS/fdNdy0;
222 fPtWgt = intPtS/intPt0;
223 fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
225 fYWgt = intYS/fdNdy0;
226 fPtWgt = (fPtMax-fPtMin)/intPt0;
227 fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
230 // particle decay related initialization
231 fDecayer->SetForceDecay(fForceDecay);
238 //____________________________________________________________
239 void AliGenParam::Generate()
242 // Generate 'npart' of light and heavy mesons (J/Psi, upsilon or phi, Pion,
243 // Kaons, Etas, Omegas) and Baryons (proton, antiprotons, neutrons and
244 // antineutrons in the the desired theta, phi and momentum windows;
245 // Gaussian smearing on the vertex is done if selected.
246 // The decay of heavy mesons is done using lujet,
247 // and the childern particle are tracked by GEANT
248 // However, light mesons are directly tracked by GEANT
249 // setting fForceDecay = nodecay (SetForceDecay(nodecay))
252 // Reinitialize decayer
253 fDecayer->SetForceDecay(fForceDecay);
257 Float_t polar[3]= {0,0,0}; // Polarisation of the parent particle (for GEANT tracking)
258 Float_t origin0[3]; // Origin of the generated parent particle (for GEANT tracking)
259 Float_t pt, pl, ptot; // Transverse, logitudinal and total momenta of the parent particle
260 Float_t phi, theta; // Phi and theta spherical angles of the parent particle momentum
262 och[3]; // Momentum, polarisation and origin of the children particles from lujet
267 static TClonesArray *particles;
269 if(!particles) particles = new TClonesArray("TParticle",1000);
271 TDatabasePDG *pDataBase = TDatabasePDG::Instance();
275 // Calculating vertex position per event
276 for (j=0;j<3;j++) origin0[j]=fOrigin[j];
277 if(fVertexSmear==kPerEvent) {
279 for (j=0;j<3;j++) origin0[j]=fVertex[j];
284 // Generating fNpart particles
289 Int_t iPart = fIpParaFunc(fRandom);
290 fChildWeight=(fDecayer->GetPartialBranchingRatio(iPart))*fParentWeight;
291 TParticlePDG *particle = pDataBase->GetParticle(iPart);
292 Float_t am = particle->Mass();
297 phi=fPhiMin+random[0]*(fPhiMax-fPhiMin);
300 ty = TMath::TanH(fYPara->GetRandom());
303 if (fAnalog == kAnalog) {
304 pt=fPtPara->GetRandom();
308 pt=fPtMin+random[1]*(fPtMax-fPtMin);
310 wgtp=fParentWeight*fPtParaFunc(& ptd, &dummy);
311 wgtch=fChildWeight*fPtParaFunc(& ptd, &dummy);
313 xmt=sqrt(pt*pt+am*am);
314 if (TMath::Abs(ty)==1.) {
317 "Division by 0: Please check you rapidity range !");
320 pl=xmt*ty/sqrt(1.-ty*ty);
321 theta=TMath::ATan2(pt,pl);
323 if(theta<fThetaMin || theta>fThetaMax) continue;
324 ptot=TMath::Sqrt(pt*pt+pl*pl);
326 if(ptot<fPMin || ptot>fPMax) continue;
328 p[0]=pt*TMath::Cos(phi);
329 p[1]=pt*TMath::Sin(phi);
331 if(fVertexSmear==kPerTrack) {
335 fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
336 TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
340 // Looking at fForceDecay :
341 // if fForceDecay != none Primary particle decays using
342 // AliPythia and children are tracked by GEANT
344 // if fForceDecay == none Primary particle is tracked by GEANT
345 // (In the latest, make sure that GEANT actually does all the decays you want)
348 if (fForceDecay != kNoDecay) {
349 // Using lujet to decay particle
350 Float_t energy=TMath::Sqrt(ptot*ptot+am*am);
351 TLorentzVector pmom(p[0], p[1], p[2], energy);
352 fDecayer->Decay(iPart,&pmom);
354 // select decay particles
355 Int_t np=fDecayer->ImportParticles(particles);
357 // Selecting GeometryAcceptance for particles fPdgCodeParticleforAcceptanceCut;
358 if (fGeometryAcceptance)
359 if (!CheckAcceptanceGeometry(np,particles)) continue;
361 Int_t* pFlag = new Int_t[np];
362 Int_t* pParent = new Int_t[np];
363 Int_t* pSelected = new Int_t[np];
364 Int_t* trackIt = new Int_t[np];
366 for (i=0; i<np; i++) {
373 TParticle* iparticle = (TParticle *) particles->At(0);
375 for (i = 1; i<np ; i++) {
377 iparticle = (TParticle *) particles->At(i);
378 Int_t kf = iparticle->GetPdgCode();
379 Int_t ks = iparticle->GetStatusCode();
383 ipF = iparticle->GetFirstDaughter();
384 ipL = iparticle->GetLastDaughter();
385 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
389 // flag decay products of particles with long life-time (c tau > .3 mum)
392 // TParticlePDG *particle = pDataBase->GetParticle(kf);
394 Double_t lifeTime = fDecayer->GetLifetime(kf);
395 // Double_t mass = particle->Mass();
396 // Double_t width = particle->Width();
397 if (lifeTime > (Double_t) fMaxLifeTime) {
398 ipF = iparticle->GetFirstDaughter();
399 ipL = iparticle->GetLastDaughter();
400 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
409 if (ChildSelected(TMath::Abs(kf)) || fForceDecay == kAll && trackIt[i])
412 pc[0]=iparticle->Px();
413 pc[1]=iparticle->Py();
414 pc[2]=iparticle->Pz();
415 Bool_t childok = KinematicSelection(iparticle, 1);
426 } // if child selection
428 } // decay particle loop
429 } // if decay products
432 if ((fCutOnChild && ncsel >0) || !fCutOnChild){
436 PushTrack(0, -1, iPart, p, origin0, polar, 0, kPPrimary, nt, wgtp);
442 for (i = 1; i < np; i++) {
444 TParticle* iparticle = (TParticle *) particles->At(i);
445 Int_t kf = iparticle->GetPdgCode();
446 Int_t ipa = iparticle->GetFirstMother()-1;
448 och[0] = origin0[0]+iparticle->Vx()/10;
449 och[1] = origin0[1]+iparticle->Vy()/10;
450 och[2] = origin0[2]+iparticle->Vz()/10;
451 pc[0] = iparticle->Px();
452 pc[1] = iparticle->Py();
453 pc[2] = iparticle->Pz();
456 iparent = pParent[ipa];
461 PushTrack(fTrackIt*trackIt[i], iparent, kf,
463 0, kPDecay, nt, wgtch);
470 if (pFlag) delete[] pFlag;
471 if (pParent) delete[] pParent;
472 if (pSelected) delete[] pSelected;
473 if (trackIt) delete[] trackIt;
474 } // kinematic selection
475 else // nodecay option, so parent will be tracked by GEANT (pions, kaons, eta, omegas, baryons)
478 PushTrack(fTrackIt,-1,iPart,p,origin0,polar,0,kPPrimary,nt,wgtp);
484 SetHighWaterMark(nt);
486 //____________________________________________________________________________________
487 Float_t AliGenParam::GetRelativeArea(Float_t ptMin, Float_t ptMax, Float_t yMin, Float_t yMax, Float_t phiMin, Float_t phiMax)
490 // Normalisation for selected kinematic region
493 fPtPara->Integral(ptMin,ptMax,(Double_t *)0,1.e-6) / fPtPara->Integral( fPtPara->GetXmin(), fPtPara->GetXmax(),(Double_t *)0,1.e-6) *
494 fYPara->Integral(yMin,yMax,(Double_t *)0,1.e-6)/fYPara->Integral(fYPara->GetXmin(),fYPara->GetXmax(),(Double_t *)0,1.e-6) *
495 (phiMax-phiMin)/360.;
496 return TMath::Abs(ratio);
499 //____________________________________________________________________________________
501 void AliGenParam::Draw( const char * /*opt*/)
504 // Draw the pT and y Distributions
506 TCanvas *c0 = new TCanvas("c0","Canvas 0",400,10,600,700);
510 fPtPara->GetHistogram()->SetXTitle("p_{T} (GeV)");
513 fYPara->GetHistogram()->SetXTitle("y");