]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenParam.cxx
No TF1/TF2 objects with the same name allowed in Root v3-10-01a
[u/mrichter/AliRoot.git] / EVGEN / AliGenParam.cxx
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 "AliGenParam.h"
25 #include "AliDecayer.h"
26 #include "AliGenMUONlib.h"
27 #include "AliRun.h"
28 #include <TParticle.h>
29 #include <TParticlePDG.h>
30 #include <TDatabasePDG.h>
31 #include <TLorentzVector.h>
32 #include <TVirtualMC.h>
33
34 #include <TF1.h>
35 #include <TCanvas.h>
36 #include <TH1.h>
37 #include "AliMC.h"
38
39 ClassImp(AliGenParam)
40
41 //------------------------------------------------------------
42
43   //Begin_Html
44   /*
45     <img src="picts/AliGenParam.gif">
46   */
47   //End_Html
48
49 //____________________________________________________________
50 //____________________________________________________________
51 AliGenParam::AliGenParam()
52 {
53 // Deafault constructor
54     fPtPara = 0;
55     fYPara  = 0;
56     fParam  = 0;
57     fAnalog = kAnalog;
58     SetDeltaPt();
59     fDecayer = 0;
60
61
62 }
63
64 AliGenParam::AliGenParam(Int_t npart, AliGenLib * Library,  Int_t param, char* tname):AliGenMC(npart)
65 {
66 // Constructor using number of particles parameterisation id and library
67     fName = "Param";
68     fTitle= "Particle Generator using pT and y parameterisation";
69     
70     fPtParaFunc = Library->GetPt(param, tname);
71     fYParaFunc  = Library->GetY (param, tname);
72     fIpParaFunc = Library->GetIp(param, tname);
73     
74     fPtPara = 0;
75     fYPara  = 0;
76     fParam  = param;
77     fAnalog = kAnalog;
78     SetForceDecay();
79     SetDeltaPt(); 
80 }
81
82 //____________________________________________________________
83
84 AliGenParam::AliGenParam(Int_t npart, Int_t param, char* tname):AliGenMC(npart)
85 {
86 // Constructor using parameterisation id and number of particles
87 //
88     fName = "Param";
89     fTitle= "Particle Generator using pT and y parameterisation";
90       
91     AliGenLib* pLibrary = new AliGenMUONlib();
92  
93     fPtParaFunc = pLibrary->GetPt(param, tname);
94     fYParaFunc  = pLibrary->GetY (param, tname);
95     fIpParaFunc = pLibrary->GetIp(param, tname);
96     
97     fPtPara = 0;
98     fYPara  = 0;
99     fParam  = param;
100     fAnalog = kAnalog;
101     fChildSelect.Set(5);
102     for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
103     SetForceDecay();
104     SetCutOnChild();
105     SetChildMomentumRange();
106     SetChildPtRange();
107     SetChildPhiRange();
108     SetChildThetaRange(); 
109     SetDeltaPt(); 
110 }
111
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 *))                 
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
143 AliGenParam::AliGenParam(const AliGenParam & Param)
144     :AliGenMC(Param)
145 {
146 // Copy constructor
147     Param.Copy(*this);
148 }
149
150 //____________________________________________________________
151 AliGenParam::~AliGenParam()
152 {
153 // Destructor
154     delete  fPtPara;
155     delete  fYPara;
156 }
157
158 //____________________________________________________________
159 void 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  
170     fPtPara = new TF1("Pt-Parametrization",fPtParaFunc,fPtMin,fPtMax,0);
171 //  Set representation precision to 10 MeV
172     Int_t npx= Int_t((fPtMax-fPtMin)/fDeltaPt);
173     
174     fPtPara->SetNpx(npx);
175     
176     fYPara  = new TF1("Y-Parametrization",fYParaFunc,fYMin,fYMax,0);
177     TF1 ptPara("Pt-Parametrization(0,15)",fPtParaFunc,0,15,0);
178     TF1 yPara("Y-Parametrization(-6,6)",fYParaFunc,-6,6,0);
179
180 //
181 // dN/dy| y=0
182     Double_t y1=0;
183     Double_t y2=0;
184     
185     fdNdy0=fYParaFunc(&y1,&y2);
186 //
187 // Integral over generation region
188     Float_t intYS  = yPara.Integral(fYMin, fYMax);
189     Float_t intPt0 = ptPara.Integral(0,15);
190     Float_t intPtS = ptPara.Integral(fPtMin,fPtMax);
191     Float_t phiWgt=(fPhiMax-fPhiMin)/2./TMath::Pi();
192     if (fAnalog == kAnalog) {
193         fYWgt  = intYS/fdNdy0;
194         fPtWgt = intPtS/intPt0;
195         fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
196     } else {
197         fYWgt = intYS/fdNdy0;
198         fPtWgt = (fPtMax-fPtMin)/intPt0;
199         fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
200     }
201 //
202 // particle decay related initialization
203     fDecayer->SetForceDecay(fForceDecay);
204     fDecayer->Init();
205
206 //
207     AliGenMC::Init();
208 }
209
210 //____________________________________________________________
211 void AliGenParam::Generate()
212 {
213 //
214 // Generate 'npart' of light and heavy mesons (J/Psi, upsilon or phi, Pion,
215 // Kaons, Etas, Omegas) and Baryons (proton, antiprotons, neutrons and 
216 // antineutrons in the the desired theta, phi and momentum windows; 
217 // Gaussian smearing on the vertex is done if selected. 
218 // The decay of heavy mesons is done using lujet, 
219 //    and the childern particle are tracked by GEANT
220 // However, light mesons are directly tracked by GEANT 
221 // setting fForceDecay = nodecay (SetForceDecay(nodecay)) 
222 //
223 //
224 //  Reinitialize decayer
225     fDecayer->Init();
226 //
227   Float_t polar[3]= {0,0,0};  // Polarisation of the parent particle (for GEANT tracking)
228   Float_t origin0[3];         // Origin of the generated parent particle (for GEANT tracking)
229   Float_t pt, pl, ptot;       // Transverse, logitudinal and total momenta of the parent particle
230   Float_t phi, theta;         // Phi and theta spherical angles of the parent particle momentum
231   Float_t p[3], pc[3], 
232           och[3];             // Momentum, polarisation and origin of the children particles from lujet
233   Double_t ty, xmt;
234   Int_t nt, i, j;
235   Float_t  wgtp, wgtch;
236   Double_t dummy;
237   static TClonesArray *particles;
238   //
239   if(!particles) particles = new TClonesArray("TParticle",1000);
240   
241   TDatabasePDG *pDataBase = TDatabasePDG::Instance();
242   //
243   Float_t random[6];
244  
245 // Calculating vertex position per event
246   for (j=0;j<3;j++) origin0[j]=fOrigin[j];
247   if(fVertexSmear==kPerEvent) {
248       Vertex();
249       for (j=0;j<3;j++) origin0[j]=fVertex[j];
250   }
251   
252   Int_t ipa=0;
253   
254 // Generating fNpart particles
255   while (ipa<fNpart) {
256       while(1) {
257 //
258 // particle type
259           Int_t iPart = fIpParaFunc(fRandom);
260           fChildWeight=(fDecayer->GetPartialBranchingRatio(iPart))*fParentWeight;          
261           TParticlePDG *particle = pDataBase->GetParticle(iPart);
262           Float_t am = particle->Mass();
263           
264           Rndm(random,2);
265 //
266 // phi
267           phi=fPhiMin+random[0]*(fPhiMax-fPhiMin);
268 //
269 // y
270           ty = TMath::TanH(fYPara->GetRandom());
271 //
272 // pT
273           if (fAnalog == kAnalog) {
274               pt=fPtPara->GetRandom();
275               wgtp=fParentWeight;
276               wgtch=fChildWeight;
277           } else {
278               pt=fPtMin+random[1]*(fPtMax-fPtMin);
279               Double_t ptd=pt;
280               wgtp=fParentWeight*fPtParaFunc(& ptd, &dummy);
281               wgtch=fChildWeight*fPtParaFunc(& ptd, &dummy);
282           }
283           xmt=sqrt(pt*pt+am*am);
284           if (TMath::Abs(ty)==1.) {
285               ty=0.;
286               Fatal("AliGenParam", 
287                     "Division by 0: Please check you rapidity range !");
288           }
289           
290           pl=xmt*ty/sqrt(1.-ty*ty);
291           theta=TMath::ATan2(pt,pl);
292 // Cut on theta
293           if(theta<fThetaMin || theta>fThetaMax) continue;
294           ptot=TMath::Sqrt(pt*pt+pl*pl);
295 // Cut on momentum
296           if(ptot<fPMin || ptot>fPMax) continue;
297 //
298           p[0]=pt*TMath::Cos(phi);
299           p[1]=pt*TMath::Sin(phi);
300           p[2]=pl;
301           if(fVertexSmear==kPerTrack) {
302               Rndm(random,6);
303               for (j=0;j<3;j++) {
304                   origin0[j]=
305                       fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
306                       TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
307               }
308           }
309           
310 // Looking at fForceDecay : 
311 // if fForceDecay != none Primary particle decays using 
312 // AliPythia and children are tracked by GEANT
313 //
314 // if fForceDecay == none Primary particle is tracked by GEANT 
315 // (In the latest, make sure that GEANT actually does all the decays you want)    
316 //
317
318           if (fForceDecay != kNoDecay) {
319 // Using lujet to decay particle
320               Float_t energy=TMath::Sqrt(ptot*ptot+am*am);
321               TLorentzVector pmom(p[0], p[1], p[2], energy);
322               fDecayer->Decay(iPart,&pmom);
323 //
324 // select decay particles
325               Int_t np=fDecayer->ImportParticles(particles);
326
327               //  Selecting  GeometryAcceptance for particles fPdgCodeParticleforAcceptanceCut;
328               if (fGeometryAcceptance) 
329                 if (!CheckAcceptanceGeometry(np,particles)) continue;
330               Int_t ncsel=0;
331               Int_t* pFlag      = new Int_t[np];
332               Int_t* pParent    = new Int_t[np];
333               Int_t* pSelected  = new Int_t[np];
334               Int_t* trackIt    = new Int_t[np];
335
336               for (i=0; i<np; i++) {
337                   pFlag[i]     =  0;
338                   pSelected[i] =  0;
339                   pParent[i]   = -1;
340               }
341               
342               if (np >1) {
343                   TParticle* iparticle =  (TParticle *) particles->At(0);
344                   Int_t ipF, ipL;
345                   for (i = 1; i<np ; i++) {
346                       trackIt[i] = 1;
347                       iparticle = (TParticle *) particles->At(i);
348                       Int_t kf = iparticle->GetPdgCode();
349                       Int_t ks = iparticle->GetStatusCode();
350 // flagged particle
351
352                       if (pFlag[i] == 1) {
353                           ipF = iparticle->GetFirstDaughter();
354                           ipL = iparticle->GetLastDaughter();   
355                           if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
356                           continue;
357                       }
358
359 // flag decay products of particles with long life-time (c tau > .3 mum)                      
360                       
361                       if (ks != 1) { 
362 //                        TParticlePDG *particle = pDataBase->GetParticle(kf);
363                           
364                           Double_t lifeTime = fDecayer->GetLifetime(kf);
365 //                        Double_t mass     = particle->Mass();
366 //                        Double_t width    = particle->Width();
367                           if (lifeTime > (Double_t) fMaxLifeTime) {
368                               ipF = iparticle->GetFirstDaughter();
369                               ipL = iparticle->GetLastDaughter();       
370                               if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
371                           } else{
372                               trackIt[i]     = 0;
373                               pSelected[i]   = 1;
374                           }
375                       } // ks==1 ?
376 //
377 // children
378                       
379                       if (ChildSelected(TMath::Abs(kf)) || fForceDecay == kAll && trackIt[i])
380                       {
381                           if (fCutOnChild) {
382                               pc[0]=iparticle->Px();
383                               pc[1]=iparticle->Py();
384                               pc[2]=iparticle->Pz();
385                               Bool_t  childok = KinematicSelection(iparticle, 1);
386                               if(childok) {
387                                   pSelected[i]  = 1;
388                                   ncsel++;
389                               } else {
390                                   ncsel=-1;
391                                   break;
392                               } // child kine cuts
393                           } else {
394                               pSelected[i]  = 1;
395                               ncsel++;
396                           } // if child selection
397                       } // select muon
398                   } // decay particle loop
399               } // if decay products
400               
401               Int_t iparent;
402               if ((fCutOnChild && ncsel >0) || !fCutOnChild){
403                   ipa++;
404 //
405 // Parent
406                   PushTrack(0, -1, iPart, p, origin0, polar, 0, kPPrimary, nt, wgtp);
407                   pParent[0] = nt;
408                   KeepTrack(nt); 
409 //
410 // Decay Products
411 //                
412                   for (i = 1; i < np; i++) {
413                       if (pSelected[i]) {
414                           TParticle* iparticle = (TParticle *) particles->At(i);
415                           Int_t kf  = iparticle->GetPdgCode();
416                           Int_t ipa = iparticle->GetFirstMother()-1;
417                           
418                           och[0] = origin0[0]+iparticle->Vx()/10;
419                           och[1] = origin0[1]+iparticle->Vy()/10;
420                           och[2] = origin0[2]+iparticle->Vz()/10;
421                           pc[0]  = iparticle->Px();
422                           pc[1]  = iparticle->Py();
423                           pc[2]  = iparticle->Pz();
424                           
425                           if (ipa > -1) {
426                               iparent = pParent[ipa];
427                           } else {
428                               iparent = -1;
429                           }
430                          
431                           PushTrack(fTrackIt*trackIt[i], iparent, kf,
432                                            pc, och, polar,
433                                            0, kPDecay, nt, wgtch);
434                           pParent[i] = nt;
435                           KeepTrack(nt); 
436                       } // Selected
437                   } // Particle loop 
438               }  // Decays by Lujet
439               particles->Clear();
440               if (pFlag)      delete[] pFlag;
441               if (pParent)    delete[] pParent;
442               if (pSelected)  delete[] pSelected;          
443               if (trackIt)    delete[] trackIt;
444           } // kinematic selection
445           else  // nodecay option, so parent will be tracked by GEANT (pions, kaons, eta, omegas, baryons)
446           {
447             gAlice->GetMCApp()->
448                 PushTrack(fTrackIt,-1,iPart,p,origin0,polar,0,kPPrimary,nt,wgtp);
449             ipa++; 
450           }
451           break;
452     } // while
453   } // event loop
454   SetHighWaterMark(nt);
455 }
456
457 void AliGenParam::Draw()
458 {
459     //
460     // Draw the pT and y Distributions
461     //
462      TCanvas *c0 = new TCanvas("c0","Canvas 0",400,10,600,700);
463      c0->Divide(2,1);
464      c0->cd(1);
465      fPtPara->Draw();
466      fPtPara->GetHistogram()->SetXTitle("p_{T} (GeV)");     
467      c0->cd(2);
468      fYPara->Draw();
469      fYPara->GetHistogram()->SetXTitle("y");     
470 }
471
472 AliGenParam& AliGenParam::operator=(const  AliGenParam& rhs)
473 {
474 // Assignment operator
475     rhs.Copy(*this);
476     return *this;
477 }
478
479
480