]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - EVGEN/AliGenParam.cxx
Corrected path to the default CDB storage (A.Colla)
[u/mrichter/AliRoot.git] / EVGEN / AliGenParam.cxx
... / ...
CommitLineData
1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
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 **************************************************************************/
15
16/* $Id$ */
17
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
23
24#include <TCanvas.h>
25#include <TDatabasePDG.h>
26#include <TF1.h>
27#include <TH1F.h>
28#include <TLorentzVector.h>
29#include <TMath.h>
30#include <TParticle.h>
31#include <TParticlePDG.h>
32#include <TROOT.h>
33#include <TVirtualMC.h>
34
35#include "AliDecayer.h"
36#include "AliGenMUONlib.h"
37#include "AliGenParam.h"
38#include "AliMC.h"
39#include "AliRun.h"
40
41ClassImp(AliGenParam)
42
43//------------------------------------------------------------
44
45 //Begin_Html
46 /*
47 <img src="picts/AliGenParam.gif">
48 */
49 //End_Html
50
51//____________________________________________________________
52AliGenParam::AliGenParam()
53{
54// Deafault constructor
55 fPtPara = 0;
56 fYPara = 0;
57 fParam = 0;
58 fAnalog = kAnalog;
59 SetDeltaPt();
60 fDecayer = 0;
61
62
63}
64//____________________________________________________________
65AliGenParam::AliGenParam(Int_t npart, AliGenLib * Library, Int_t param, char* tname):AliGenMC(npart)
66{
67// Constructor using number of particles parameterisation id and library
68 fName = "Param";
69 fTitle= "Particle Generator using pT and y parameterisation";
70
71 fPtParaFunc = Library->GetPt(param, tname);
72 fYParaFunc = Library->GetY (param, tname);
73 fIpParaFunc = Library->GetIp(param, tname);
74
75 fPtPara = 0;
76 fYPara = 0;
77 fParam = param;
78 fAnalog = kAnalog;
79 SetForceDecay();
80 SetDeltaPt();
81}
82//____________________________________________________________
83AliGenParam::AliGenParam(Int_t npart, Int_t param, const char* tname, const char* name):AliGenMC(npart)
84{
85// Constructor using parameterisation id and number of particles
86//
87 fName = name;
88 fTitle= "Particle Generator using pT and y parameterisation";
89
90 AliGenLib* pLibrary = new AliGenMUONlib();
91
92 fPtParaFunc = pLibrary->GetPt(param, tname);
93 fYParaFunc = pLibrary->GetY (param, tname);
94 fIpParaFunc = pLibrary->GetIp(param, tname);
95
96 fPtPara = 0;
97 fYPara = 0;
98 fParam = param;
99 fAnalog = kAnalog;
100 fChildSelect.Set(5);
101 for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
102 SetForceDecay();
103 SetCutOnChild();
104 SetChildMomentumRange();
105 SetChildPtRange();
106 SetChildPhiRange();
107 SetChildThetaRange();
108 SetDeltaPt();
109}
110//____________________________________________________________
111
112AliGenParam::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 *))
116 :AliGenMC(npart)
117{
118// Constructor
119// Gines Martinez 1/10/99
120 fName = "Param";
121 fTitle= "Particle Generator using pT and y parameterisation";
122
123 fPtParaFunc = PtPara;
124 fYParaFunc = YPara;
125 fIpParaFunc = IpPara;
126//
127 fPtPara = 0;
128 fYPara = 0;
129 fParam = param;
130 fAnalog = kAnalog;
131 fChildSelect.Set(5);
132 for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
133 SetForceDecay();
134 SetCutOnChild();
135 SetChildMomentumRange();
136 SetChildPtRange();
137 SetChildPhiRange();
138 SetChildThetaRange();
139 SetDeltaPt();
140}
141
142
143AliGenParam::AliGenParam(const AliGenParam & Param)
144 :AliGenMC(Param)
145{
146// Copy constructor
147 Param.Copy(*this);
148}
149
150//____________________________________________________________
151AliGenParam::~AliGenParam()
152{
153// Destructor
154 delete fPtPara;
155 delete fYPara;
156}
157
158//____________________________________________________________
159void AliGenParam::Init()
160{
161// Initialisation
162
163 if (gMC) fDecayer = gMC->GetDecayer();
164 //Begin_Html
165 /*
166 <img src="picts/AliGenParam.gif">
167 */
168 //End_Html
169 char name[256];
170 sprintf(name, "pt-parameterisation for %s", GetName());
171
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);
177
178 fPtPara->SetNpx(npx);
179
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);
184
185
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);
190
191//
192// dN/dy| y=0
193 Double_t y1=0;
194 Double_t y2=0;
195
196 fdNdy0=fYParaFunc(&y1,&y2);
197//
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;
207 } else {
208 fYWgt = intYS/fdNdy0;
209 fPtWgt = (fPtMax-fPtMin)/intPt0;
210 fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
211 }
212//
213// particle decay related initialization
214 fDecayer->SetForceDecay(fForceDecay);
215 fDecayer->Init();
216
217//
218 AliGenMC::Init();
219}
220
221//____________________________________________________________
222void AliGenParam::Generate()
223{
224//
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))
233//
234//
235// Reinitialize decayer
236 fDecayer->SetForceDecay(fForceDecay);
237 fDecayer->Init();
238
239//
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
244 Float_t p[3], pc[3],
245 och[3]; // Momentum, polarisation and origin of the children particles from lujet
246 Double_t ty, xmt;
247 Int_t nt, i, j;
248 Float_t wgtp, wgtch;
249 Double_t dummy;
250 static TClonesArray *particles;
251 //
252 if(!particles) particles = new TClonesArray("TParticle",1000);
253
254 TDatabasePDG *pDataBase = TDatabasePDG::Instance();
255 //
256 Float_t random[6];
257
258// Calculating vertex position per event
259 for (j=0;j<3;j++) origin0[j]=fOrigin[j];
260 if(fVertexSmear==kPerEvent) {
261 Vertex();
262 for (j=0;j<3;j++) origin0[j]=fVertex[j];
263 }
264
265 Int_t ipa=0;
266
267// Generating fNpart particles
268 while (ipa<fNpart) {
269 while(1) {
270//
271// particle type
272 Int_t iPart = fIpParaFunc(fRandom);
273 fChildWeight=(fDecayer->GetPartialBranchingRatio(iPart))*fParentWeight;
274 TParticlePDG *particle = pDataBase->GetParticle(iPart);
275 Float_t am = particle->Mass();
276
277 Rndm(random,2);
278//
279// phi
280 phi=fPhiMin+random[0]*(fPhiMax-fPhiMin);
281//
282// y
283 ty = TMath::TanH(fYPara->GetRandom());
284//
285// pT
286 if (fAnalog == kAnalog) {
287 pt=fPtPara->GetRandom();
288 wgtp=fParentWeight;
289 wgtch=fChildWeight;
290 } else {
291 pt=fPtMin+random[1]*(fPtMax-fPtMin);
292 Double_t ptd=pt;
293 wgtp=fParentWeight*fPtParaFunc(& ptd, &dummy);
294 wgtch=fChildWeight*fPtParaFunc(& ptd, &dummy);
295 }
296 xmt=sqrt(pt*pt+am*am);
297 if (TMath::Abs(ty)==1.) {
298 ty=0.;
299 Fatal("AliGenParam",
300 "Division by 0: Please check you rapidity range !");
301 }
302
303 pl=xmt*ty/sqrt(1.-ty*ty);
304 theta=TMath::ATan2(pt,pl);
305// Cut on theta
306 if(theta<fThetaMin || theta>fThetaMax) continue;
307 ptot=TMath::Sqrt(pt*pt+pl*pl);
308// Cut on momentum
309 if(ptot<fPMin || ptot>fPMax) continue;
310//
311 p[0]=pt*TMath::Cos(phi);
312 p[1]=pt*TMath::Sin(phi);
313 p[2]=pl;
314 if(fVertexSmear==kPerTrack) {
315 Rndm(random,6);
316 for (j=0;j<3;j++) {
317 origin0[j]=
318 fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
319 TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
320 }
321 }
322
323// Looking at fForceDecay :
324// if fForceDecay != none Primary particle decays using
325// AliPythia and children are tracked by GEANT
326//
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)
329//
330
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);
336//
337// select decay particles
338 Int_t np=fDecayer->ImportParticles(particles);
339
340 // Selecting GeometryAcceptance for particles fPdgCodeParticleforAcceptanceCut;
341 if (fGeometryAcceptance)
342 if (!CheckAcceptanceGeometry(np,particles)) continue;
343 Int_t ncsel=0;
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];
348
349 for (i=0; i<np; i++) {
350 pFlag[i] = 0;
351 pSelected[i] = 0;
352 pParent[i] = -1;
353 }
354
355 if (np >1) {
356 TParticle* iparticle = (TParticle *) particles->At(0);
357 Int_t ipF, ipL;
358 for (i = 1; i<np ; i++) {
359 trackIt[i] = 1;
360 iparticle = (TParticle *) particles->At(i);
361 Int_t kf = iparticle->GetPdgCode();
362 Int_t ks = iparticle->GetStatusCode();
363// flagged particle
364
365 if (pFlag[i] == 1) {
366 ipF = iparticle->GetFirstDaughter();
367 ipL = iparticle->GetLastDaughter();
368 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
369 continue;
370 }
371
372// flag decay products of particles with long life-time (c tau > .3 mum)
373
374 if (ks != 1) {
375// TParticlePDG *particle = pDataBase->GetParticle(kf);
376
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;
384 } else{
385 trackIt[i] = 0;
386 pSelected[i] = 1;
387 }
388 } // ks==1 ?
389//
390// children
391
392 if (ChildSelected(TMath::Abs(kf)) || fForceDecay == kAll && trackIt[i])
393 {
394 if (fCutOnChild) {
395 pc[0]=iparticle->Px();
396 pc[1]=iparticle->Py();
397 pc[2]=iparticle->Pz();
398 Bool_t childok = KinematicSelection(iparticle, 1);
399 if(childok) {
400 pSelected[i] = 1;
401 ncsel++;
402 } else {
403 ncsel=-1;
404 break;
405 } // child kine cuts
406 } else {
407 pSelected[i] = 1;
408 ncsel++;
409 } // if child selection
410 } // select muon
411 } // decay particle loop
412 } // if decay products
413
414 Int_t iparent;
415 if ((fCutOnChild && ncsel >0) || !fCutOnChild){
416 ipa++;
417//
418// Parent
419 PushTrack(0, -1, iPart, p, origin0, polar, 0, kPPrimary, nt, wgtp);
420 pParent[0] = nt;
421 KeepTrack(nt);
422//
423// Decay Products
424//
425 for (i = 1; i < np; i++) {
426 if (pSelected[i]) {
427 TParticle* iparticle = (TParticle *) particles->At(i);
428 Int_t kf = iparticle->GetPdgCode();
429 Int_t ipa = iparticle->GetFirstMother()-1;
430
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();
437
438 if (ipa > -1) {
439 iparent = pParent[ipa];
440 } else {
441 iparent = -1;
442 }
443
444 PushTrack(fTrackIt*trackIt[i], iparent, kf,
445 pc, och, polar,
446 0, kPDecay, nt, wgtch);
447 pParent[i] = nt;
448 KeepTrack(nt);
449 } // Selected
450 } // Particle loop
451 } // Decays by Lujet
452 particles->Clear();
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)
459 {
460 gAlice->GetMCApp()->
461 PushTrack(fTrackIt,-1,iPart,p,origin0,polar,0,kPPrimary,nt,wgtp);
462 ipa++;
463 }
464 break;
465 } // while
466 } // event loop
467 SetHighWaterMark(nt);
468}
469//____________________________________________________________________________________
470Float_t AliGenParam::GetRelativeArea(Float_t ptMin, Float_t ptMax, Float_t yMin, Float_t yMax, Float_t phiMin, Float_t phiMax)
471{
472//
473// Normalisation for selected kinematic region
474//
475 Float_t ratio =
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);
480}
481
482//____________________________________________________________________________________
483
484void AliGenParam::Draw( const char * /*opt*/)
485{
486 //
487 // Draw the pT and y Distributions
488 //
489 TCanvas *c0 = new TCanvas("c0","Canvas 0",400,10,600,700);
490 c0->Divide(2,1);
491 c0->cd(1);
492 fPtPara->Draw();
493 fPtPara->GetHistogram()->SetXTitle("p_{T} (GeV)");
494 c0->cd(2);
495 fYPara->Draw();
496 fYPara->GetHistogram()->SetXTitle("y");
497}
498
499AliGenParam& AliGenParam::operator=(const AliGenParam& rhs)
500{
501// Assignment operator
502 rhs.Copy(*this);
503 return *this;
504}
505
506
507