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