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