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