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