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():
69 // Default constructor
71 //____________________________________________________________
72 AliGenParam::AliGenParam(Int_t npart, AliGenLib * Library, Int_t param, char* tname)
74 fPtParaFunc(Library->GetPt(param, tname)),
75 fYParaFunc (Library->GetY (param, tname)),
76 fIpParaFunc(Library->GetIp(param, tname)),
88 // Constructor using number of particles parameterisation id and library
90 fTitle= "Particle Generator using pT and y parameterisation";
94 //____________________________________________________________
95 AliGenParam::AliGenParam(Int_t npart, Int_t param, const char* tname, const char* name):
111 // Constructor using parameterisation id and number of particles
114 fTitle= "Particle Generator using pT and y parameterisation";
116 AliGenLib* pLibrary = new AliGenMUONlib();
117 fPtParaFunc = pLibrary->GetPt(param, tname);
118 fYParaFunc = pLibrary->GetY (param, tname);
119 fIpParaFunc = pLibrary->GetIp(param, tname);
123 for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
126 SetChildMomentumRange();
129 SetChildThetaRange();
131 //____________________________________________________________
133 AliGenParam::AliGenParam(Int_t npart, Int_t param,
134 Double_t (*PtPara) (Double_t*, Double_t*),
135 Double_t (*YPara ) (Double_t* ,Double_t*),
136 Int_t (*IpPara) (TRandom *))
154 // Gines Martinez 1/10/99
156 fTitle= "Particle Generator using pT and y parameterisation";
160 for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
163 SetChildMomentumRange();
166 SetChildThetaRange();
169 //____________________________________________________________
170 AliGenParam::~AliGenParam()
177 //____________________________________________________________
178 void AliGenParam::Init()
182 if (gMC) fDecayer = gMC->GetDecayer();
185 <img src="picts/AliGenParam.gif">
189 sprintf(name, "pt-parameterisation for %s", GetName());
191 if (fPtPara) fPtPara->Delete();
192 fPtPara = new TF1(name, fPtParaFunc, fPtMin, fPtMax,0);
193 gROOT->GetListOfFunctions()->Remove(fPtPara);
194 // Set representation precision to 10 MeV
195 Int_t npx= Int_t((fPtMax - fPtMin) / fDeltaPt);
197 fPtPara->SetNpx(npx);
199 sprintf(name, "y-parameterisation for %s", GetName());
200 if (fYPara) fYPara->Delete();
201 fYPara = new TF1(name, fYParaFunc, fYMin, fYMax, 0);
202 gROOT->GetListOfFunctions()->Remove(fYPara);
205 sprintf(name, "pt-for-%s", GetName());
206 TF1 ptPara(name ,fPtParaFunc, 0, 15, 0);
207 sprintf(name, "y-for-%s", GetName());
208 TF1 yPara(name, fYParaFunc, -6, 6, 0);
215 fdNdy0=fYParaFunc(&y1,&y2);
217 // Integral over generation region
218 Float_t intYS = yPara.Integral(fYMin, fYMax,(Double_t*) 0x0,1.e-6);
219 Float_t intPt0 = ptPara.Integral(0,15,(Double_t *) 0x0,1.e-6);
220 Float_t intPtS = ptPara.Integral(fPtMin,fPtMax,(Double_t*) 0x0,1.e-6);
221 Float_t phiWgt=(fPhiMax-fPhiMin)/2./TMath::Pi();
222 if (fAnalog == kAnalog) {
223 fYWgt = intYS/fdNdy0;
224 fPtWgt = intPtS/intPt0;
225 fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
227 fYWgt = intYS/fdNdy0;
228 fPtWgt = (fPtMax-fPtMin)/intPt0;
229 fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
232 // particle decay related initialization
233 fDecayer->SetForceDecay(fForceDecay);
240 //____________________________________________________________
241 void AliGenParam::Generate()
244 // Generate 'npart' of light and heavy mesons (J/Psi, upsilon or phi, Pion,
245 // Kaons, Etas, Omegas) and Baryons (proton, antiprotons, neutrons and
246 // antineutrons in the the desired theta, phi and momentum windows;
247 // Gaussian smearing on the vertex is done if selected.
248 // The decay of heavy mesons is done using lujet,
249 // and the childern particle are tracked by GEANT
250 // However, light mesons are directly tracked by GEANT
251 // setting fForceDecay = nodecay (SetForceDecay(nodecay))
254 // Reinitialize decayer
255 fDecayer->SetForceDecay(fForceDecay);
259 Float_t polar[3]= {0,0,0}; // Polarisation of the parent particle (for GEANT tracking)
260 Float_t origin0[3]; // Origin of the generated parent particle (for GEANT tracking)
261 Float_t pt, pl, ptot; // Transverse, logitudinal and total momenta of the parent particle
262 Float_t phi, theta; // Phi and theta spherical angles of the parent particle momentum
264 och[3]; // Momentum, polarisation and origin of the children particles from lujet
269 static TClonesArray *particles;
271 if(!particles) particles = new TClonesArray("TParticle",1000);
273 TDatabasePDG *pDataBase = TDatabasePDG::Instance();
277 // Calculating vertex position per event
278 for (j=0;j<3;j++) origin0[j]=fOrigin[j];
279 if(fVertexSmear==kPerEvent) {
281 for (j=0;j<3;j++) origin0[j]=fVertex[j];
286 // Generating fNpart particles
293 Int_t iPart = fIpParaFunc(fRandom);
294 fChildWeight=(fDecayer->GetPartialBranchingRatio(iPart))*fParentWeight;
295 TParticlePDG *particle = pDataBase->GetParticle(iPart);
296 Float_t am = particle->Mass();
301 phi=fPhiMin+random[0]*(fPhiMax-fPhiMin);
304 ty = TMath::TanH(fYPara->GetRandom());
307 if (fAnalog == kAnalog) {
308 pt=fPtPara->GetRandom();
312 pt=fPtMin+random[1]*(fPtMax-fPtMin);
314 wgtp=fParentWeight*fPtParaFunc(& ptd, &dummy);
315 wgtch=fChildWeight*fPtParaFunc(& ptd, &dummy);
317 xmt=sqrt(pt*pt+am*am);
318 if (TMath::Abs(ty)==1.) {
321 "Division by 0: Please check you rapidity range !");
324 pl=xmt*ty/sqrt(1.-ty*ty);
325 theta=TMath::ATan2(pt,pl);
327 if(theta<fThetaMin || theta>fThetaMax) continue;
328 ptot=TMath::Sqrt(pt*pt+pl*pl);
330 if(ptot<fPMin || ptot>fPMax) continue;
332 p[0]=pt*TMath::Cos(phi);
333 p[1]=pt*TMath::Sin(phi);
335 if(fVertexSmear==kPerTrack) {
339 fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
340 TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
344 // Looking at fForceDecay :
345 // if fForceDecay != none Primary particle decays using
346 // AliPythia and children are tracked by GEANT
348 // if fForceDecay == none Primary particle is tracked by GEANT
349 // (In the latest, make sure that GEANT actually does all the decays you want)
352 if (fForceDecay != kNoDecay) {
353 // Using lujet to decay particle
354 Float_t energy=TMath::Sqrt(ptot*ptot+am*am);
355 TLorentzVector pmom(p[0], p[1], p[2], energy);
356 fDecayer->Decay(iPart,&pmom);
358 // select decay particles
359 Int_t np=fDecayer->ImportParticles(particles);
361 // Selecting GeometryAcceptance for particles fPdgCodeParticleforAcceptanceCut;
362 if (fGeometryAcceptance)
363 if (!CheckAcceptanceGeometry(np,particles)) continue;
365 Int_t* pFlag = new Int_t[np];
366 Int_t* pParent = new Int_t[np];
367 Int_t* pSelected = new Int_t[np];
368 Int_t* trackIt = new Int_t[np];
370 for (i=0; i<np; i++) {
377 TParticle* iparticle = (TParticle *) particles->At(0);
379 for (i = 1; i<np ; i++) {
381 iparticle = (TParticle *) particles->At(i);
382 Int_t kf = iparticle->GetPdgCode();
383 Int_t ks = iparticle->GetStatusCode();
387 ipF = iparticle->GetFirstDaughter();
388 ipL = iparticle->GetLastDaughter();
389 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
393 // flag decay products of particles with long life-time (c tau > .3 mum)
396 // TParticlePDG *particle = pDataBase->GetParticle(kf);
398 Double_t lifeTime = fDecayer->GetLifetime(kf);
399 // Double_t mass = particle->Mass();
400 // Double_t width = particle->Width();
401 if (lifeTime > (Double_t) fMaxLifeTime) {
402 ipF = iparticle->GetFirstDaughter();
403 ipL = iparticle->GetLastDaughter();
404 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
413 if (ChildSelected(TMath::Abs(kf)) || fForceDecay == kAll && trackIt[i])
416 pc[0]=iparticle->Px();
417 pc[1]=iparticle->Py();
418 pc[2]=iparticle->Pz();
419 Bool_t childok = KinematicSelection(iparticle, 1);
430 } // if child selection
432 } // decay particle loop
433 } // if decay products
436 if ((fCutOnChild && ncsel >0) || !fCutOnChild){
440 PushTrack(0, -1, iPart, p, origin0, polar, 0, kPPrimary, nt, wgtp);
448 for (i = 1; i < np; i++) {
450 TParticle* iparticle = (TParticle *) particles->At(i);
451 Int_t kf = iparticle->GetPdgCode();
452 Int_t jpa = iparticle->GetFirstMother()-1;
454 och[0] = origin0[0]+iparticle->Vx()/10;
455 och[1] = origin0[1]+iparticle->Vy()/10;
456 och[2] = origin0[2]+iparticle->Vz()/10;
457 pc[0] = iparticle->Px();
458 pc[1] = iparticle->Py();
459 pc[2] = iparticle->Pz();
462 iparent = pParent[jpa];
467 PushTrack(fTrackIt * trackIt[i], iparent, kf,
469 0, kPDecay, nt, wgtch);
477 if (pFlag) delete[] pFlag;
478 if (pParent) delete[] pParent;
479 if (pSelected) delete[] pSelected;
480 if (trackIt) delete[] trackIt;
481 } // kinematic selection
482 else // nodecay option, so parent will be tracked by GEANT (pions, kaons, eta, omegas, baryons)
485 PushTrack(fTrackIt,-1,iPart,p,origin0,polar,0,kPPrimary,nt,wgtp);
493 SetHighWaterMark(nt);
495 AliGenEventHeader* header = new AliGenEventHeader("PARAM");
496 header->SetPrimaryVertex(fVertex);
497 header->SetNProduced(fNprimaries);
500 //____________________________________________________________________________________
501 Float_t AliGenParam::GetRelativeArea(Float_t ptMin, Float_t ptMax, Float_t yMin, Float_t yMax, Float_t phiMin, Float_t phiMax)
504 // Normalisation for selected kinematic region
507 fPtPara->Integral(ptMin,ptMax,(Double_t *)0,1.e-6) / fPtPara->Integral( fPtPara->GetXmin(), fPtPara->GetXmax(),(Double_t *)0,1.e-6) *
508 fYPara->Integral(yMin,yMax,(Double_t *)0,1.e-6)/fYPara->Integral(fYPara->GetXmin(),fYPara->GetXmax(),(Double_t *)0,1.e-6) *
509 (phiMax-phiMin)/360.;
510 return TMath::Abs(ratio);
513 //____________________________________________________________________________________
515 void AliGenParam::Draw( const char * /*opt*/)
518 // Draw the pT and y Distributions
520 TCanvas *c0 = new TCanvas("c0","Canvas 0",400,10,600,700);
524 fPtPara->GetHistogram()->SetXTitle("p_{T} (GeV)");
527 fYPara->GetHistogram()->SetXTitle("y");