]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenParam.cxx
- Possibility for vertex distribution truncation.
[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 /*
17 $Log$
18 Revision 1.31  2001/07/13 10:58:54  morsch
19 - Some coded moved to AliGenMC
20 - Improved handling of secondary vertices.
21
22 Revision 1.30  2001/06/15 07:55:04  morsch
23 Put only first generation decay products on the stack.
24
25 Revision 1.29  2001/03/27 10:58:41  morsch
26 Initialize decayer before generation. Important if run inside cocktail.
27
28 Revision 1.28  2001/03/09 13:01:41  morsch
29 - enum constants for paramterisation type (particle family) moved to AliGen*lib.h
30 - use AliGenGSIlib::kUpsilon, AliGenPHOSlib::kEtaPrime to access the constants
31
32 Revision 1.27  2001/02/02 15:21:10  morsch
33 Set high water mark after last particle.
34 Use Vertex() method for Vertex.
35
36 Revision 1.26  2000/12/21 16:24:06  morsch
37 Coding convention clean-up
38
39 Revision 1.25  2000/11/30 07:12:50  alibrary
40 Introducing new Rndm and QA classes
41
42 Revision 1.24  2000/10/18 19:11:27  hristov
43 Division by zero fixed
44
45 Revision 1.23  2000/10/02 21:28:06  fca
46 Removal of useless dependecies via forward declarations
47
48 Revision 1.22  2000/09/12 14:14:55  morsch
49 Call fDecayer->ForceDecay() at the beginning of Generate().
50
51 Revision 1.21  2000/09/08 15:39:01  morsch
52 Handle the case fForceDecay=all during the generation, i.e. select all secondaries.
53
54 Revision 1.20  2000/09/06 14:35:44  morsch
55 Use AliDecayerPythia for particle decays.
56
57 Revision 1.19  2000/07/11 18:24:56  fca
58 Coding convention corrections + few minor bug fixes
59
60 Revision 1.18  2000/06/29 21:08:27  morsch
61 All paramatrisation libraries derive from the pure virtual base class AliGenLib.
62 This allows to pass a pointer to a library directly to AliGenParam and avoids the
63 use of function pointers in Config.C.
64
65 Revision 1.17  2000/06/09 20:33:30  morsch
66 All coding rule violations except RS3 corrected
67
68 Revision 1.16  2000/05/02 07:51:31  morsch
69 - Control precision of pT sampling TF1::SetNpx(..)
70 - Correct initialisation of child-cuts in all constructors.
71 - Most coding rule violations corrected.
72
73 Revision 1.15  2000/04/03 15:42:12  morsch
74 Cuts on primary particles are separated from those on the decay products. Methods
75 SetChildMomentumRange, SetChildPtRange, SetChildPhiRange, SetChildThetaRange added.
76
77 Revision 1.14  1999/11/09 07:38:48  fca
78 Changes for compatibility with version 2.23 of ROOT
79
80 Revision 1.13  1999/11/04 11:30:31  fca
81 Correct the logics for SetForceDecay
82
83 Revision 1.12  1999/11/03 17:43:20  fca
84 New version from G.Martinez & A.Morsch
85
86 Revision 1.11  1999/09/29 09:24:14  fca
87 Introduction of the Copyright and cvs Log
88
89 */
90
91 #include "AliGenParam.h"
92 #include "AliDecayerPythia.h"
93 #include "AliGenMUONlib.h"
94 #include "AliRun.h"
95 #include <TParticle.h>
96 #include <TParticlePDG.h>
97 #include <TDatabasePDG.h>
98 #include <TLorentzVector.h>
99
100 #include <TF1.h>
101
102 ClassImp(AliGenParam)
103
104 //------------------------------------------------------------
105
106   //Begin_Html
107   /*
108     <img src="picts/AliGenParam.gif">
109   */
110   //End_Html
111
112 //____________________________________________________________
113 //____________________________________________________________
114 AliGenParam::AliGenParam()
115 {
116 // Deafault constructor
117     fPtPara = 0;
118     fYPara  = 0;
119     fParam  = 0;
120     fAnalog = kAnalog;
121     SetDeltaPt();
122 //
123 // Set random number generator   
124     sRandom = fRandom;
125 }
126
127 AliGenParam::AliGenParam(Int_t npart, AliGenLib * Library,  Int_t param, char* tname):AliGenMC(npart)
128 {
129 // Constructor using number of particles parameterisation id and library
130     
131     fPtParaFunc = Library->GetPt(param, tname);
132     fYParaFunc  = Library->GetY (param, tname);
133     fIpParaFunc = Library->GetIp(param, tname);
134     
135     fPtPara = 0;
136     fYPara  = 0;
137     fParam  = param;
138     fAnalog = kAnalog;
139     SetForceDecay();
140     SetDeltaPt(); 
141 //
142 // Set random number generator   
143     sRandom = fRandom;
144 }
145
146 //____________________________________________________________
147
148 AliGenParam::AliGenParam(Int_t npart, Int_t param, char* tname):AliGenMC(npart)
149 {
150 // Constructor using parameterisation id and number of particles
151 //      
152     AliGenLib* pLibrary = new AliGenMUONlib();
153  
154     fPtParaFunc = pLibrary->GetPt(param, tname);
155     fYParaFunc  = pLibrary->GetY (param, tname);
156     fIpParaFunc = pLibrary->GetIp(param, tname);
157     
158     fPtPara = 0;
159     fYPara  = 0;
160     fParam  = param;
161     fAnalog = kAnalog;
162     fChildSelect.Set(5);
163     for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
164     SetForceDecay();
165     SetCutOnChild();
166     SetChildMomentumRange();
167     SetChildPtRange();
168     SetChildPhiRange();
169     SetChildThetaRange(); 
170     SetDeltaPt(); 
171 }
172
173 AliGenParam::AliGenParam(Int_t npart, Int_t param,
174                          Double_t (*PtPara) (Double_t*, Double_t*),
175                          Double_t (*YPara ) (Double_t* ,Double_t*),
176                          Int_t    (*IpPara) (TRandom *))                 
177     :AliGenMC(npart)
178 {
179 // Constructor
180 // Gines Martinez 1/10/99 
181     fPtParaFunc = PtPara; 
182     fYParaFunc  = YPara;  
183     fIpParaFunc = IpPara;
184 //  
185     fPtPara = 0;
186     fYPara  = 0;
187     fParam  = param;
188     fAnalog = kAnalog;
189     fChildSelect.Set(5);
190     for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
191     SetForceDecay();
192     SetCutOnChild();
193     SetChildMomentumRange();
194     SetChildPtRange();
195     SetChildPhiRange();
196     SetChildThetaRange();  
197     SetDeltaPt();
198 }
199
200
201 AliGenParam::AliGenParam(const AliGenParam & Paramd)
202 {
203 // copy constructor
204 }
205
206 //____________________________________________________________
207 AliGenParam::~AliGenParam()
208 {
209 // Destructor
210     delete  fPtPara;
211     delete  fYPara;
212 }
213
214 //____________________________________________________________
215 void AliGenParam::Init()
216 {
217 // Initialisation
218
219     fDecayer = new AliDecayerPythia();
220   //Begin_Html
221   /*
222     <img src="picts/AliGenParam.gif">
223   */
224   //End_Html
225  
226     fPtPara = new TF1("Pt-Parametrization",fPtParaFunc,fPtMin,fPtMax,0);
227 //  Set representation precision to 10 MeV
228     Int_t npx= Int_t((fPtMax-fPtMin)/fDeltaPt);
229     
230     fPtPara->SetNpx(npx);
231     
232     fYPara  = new TF1("Y -Parametrization",fYParaFunc,fYMin,fYMax,0);
233     TF1* ptPara = new TF1("Pt-Parametrization",fPtParaFunc,0,15,0);
234     TF1* yPara  = new TF1("Y -Parametrization",fYParaFunc,-6,6,0);
235
236 //
237 // dN/dy| y=0
238     Double_t y1=0;
239     Double_t y2=0;
240     
241     fdNdy0=fYParaFunc(&y1,&y2);
242 //
243 // Integral over generation region
244     Float_t intYS  = yPara ->Integral(fYMin, fYMax);
245     Float_t intPt0 = ptPara->Integral(0,15);
246     Float_t intPtS = ptPara->Integral(fPtMin,fPtMax);
247     Float_t phiWgt=(fPhiMax-fPhiMin)/2./TMath::Pi();
248     if (fAnalog == kAnalog) {
249         fYWgt  = intYS/fdNdy0;
250         fPtWgt = intPtS/intPt0;
251         fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
252     } else {
253         fYWgt = intYS/fdNdy0;
254         fPtWgt = (fPtMax-fPtMin)/intPt0;
255         fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
256     }
257 //
258 // particle decay related initialization
259     fDecayer->SetForceDecay(fForceDecay);
260     fDecayer->Init();
261
262 //
263     AliGenMC::Init();
264 }
265
266 //____________________________________________________________
267 void AliGenParam::Generate()
268 {
269 //
270 // Generate 'npart' of light and heavy mesons (J/Psi, upsilon or phi, Pion,
271 // Kaons, Etas, Omegas) and Baryons (proton, antiprotons, neutrons and 
272 // antineutrons in the the desired theta, phi and momentum windows; 
273 // Gaussian smearing on the vertex is done if selected. 
274 // The decay of heavy mesons is done using lujet, 
275 //    and the childern particle are tracked by GEANT
276 // However, light mesons are directly tracked by GEANT 
277 // setting fForceDecay = nodecay (SetForceDecay(nodecay)) 
278 //
279 //
280 //  Reinitialize decayer
281     fDecayer->Init();
282 //
283   Float_t polar[3]= {0,0,0};  // Polarisation of the parent particle (for GEANT tracking)
284   Float_t origin0[3];         // Origin of the generated parent particle (for GEANT tracking)
285   Float_t pt, pl, ptot;       // Transverse, logitudinal and total momenta of the parent particle
286   Float_t phi, theta;         // Phi and theta spherical angles of the parent particle momentum
287   Float_t p[3], pc[3], 
288           och[3];             // Momentum, polarisation and origin of the children particles from lujet
289   Float_t ty, xmt;
290   Int_t nt, i, j;
291   Float_t  wgtp, wgtch;
292   Double_t dummy;
293   static TClonesArray *particles;
294   //
295   if(!particles) particles = new TClonesArray("TParticle",1000);
296
297   static TDatabasePDG *pDataBase = new TDatabasePDG();
298   if(!pDataBase) pDataBase = new TDatabasePDG();
299   //
300   Float_t random[6];
301  
302 // Calculating vertex position per event
303   for (j=0;j<3;j++) origin0[j]=fOrigin[j];
304   if(fVertexSmear==kPerEvent) {
305       Vertex();
306       for (j=0;j<3;j++) origin0[j]=fVertex[j];
307   }
308   
309   Int_t ipa=0;
310   
311 // Generating fNpart particles
312   while (ipa<fNpart) {
313       while(1) {
314 //
315 // particle type
316           Int_t iPart = fIpParaFunc(fRandom);
317           fChildWeight=(fDecayer->GetPartialBranchingRatio(iPart))*fParentWeight;          
318           TParticlePDG *particle = pDataBase->GetParticle(iPart);
319           Float_t am = particle->Mass();
320
321           Rndm(random,2);
322 //
323 // phi
324           phi=fPhiMin+random[0]*(fPhiMax-fPhiMin);
325 //
326 // y
327           ty=Float_t(TMath::TanH(fYPara->GetRandom()));
328 //
329 // pT
330           if (fAnalog == kAnalog) {
331               pt=fPtPara->GetRandom();
332               wgtp=fParentWeight;
333               wgtch=fChildWeight;
334           } else {
335               pt=fPtMin+random[1]*(fPtMax-fPtMin);
336               Double_t ptd=pt;
337               wgtp=fParentWeight*fPtParaFunc(& ptd, &dummy);
338               wgtch=fChildWeight*fPtParaFunc(& ptd, &dummy);
339           }
340           xmt=sqrt(pt*pt+am*am);
341       if (TMath::Abs(ty)==1) ty=0;
342           pl=xmt*ty/sqrt(1.-ty*ty);
343           theta=TMath::ATan2(pt,pl);
344 // Cut on theta
345           if(theta<fThetaMin || theta>fThetaMax) continue;
346           ptot=TMath::Sqrt(pt*pt+pl*pl);
347 // Cut on momentum
348           if(ptot<fPMin || ptot>fPMax) continue;
349 //
350           p[0]=pt*TMath::Cos(phi);
351           p[1]=pt*TMath::Sin(phi);
352           p[2]=pl;
353           if(fVertexSmear==kPerTrack) {
354               Rndm(random,6);
355               for (j=0;j<3;j++) {
356                   origin0[j]=
357                       fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
358                       TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
359               }
360           }
361           
362 // Looking at fForceDecay : 
363 // if fForceDecay != none Primary particle decays using 
364 // AliPythia and children are tracked by GEANT
365 //
366 // if fForceDecay == none Primary particle is tracked by GEANT 
367 // (In the latest, make sure that GEANT actually does all the decays you want)    
368 //
369
370           if (fForceDecay != kNoDecay) {
371 // Using lujet to decay particle
372               Float_t energy=TMath::Sqrt(ptot*ptot+am*am);
373               TLorentzVector pmom(p[0], p[1], p[2], energy);
374               fDecayer->Decay(iPart,&pmom);
375 //
376 // select decay particles
377               Int_t np=fDecayer->ImportParticles(particles);
378               Int_t ncsel=0;
379               Int_t* pFlag      = new Int_t[np];
380               Int_t* pParent    = new Int_t[np];
381               Int_t* pSelected  = new Int_t[np];
382               Int_t* trackIt    = new Int_t[np];
383
384               for (i=0; i<np; i++) {
385                   pFlag[i]     =  0;
386                   pSelected[i] =  0;
387                   pParent[i]   = -1;
388               }
389               
390               if (np >1) {
391                   TParticle* iparticle =  (TParticle *) particles->At(0);
392                   Int_t ipF, ipL;
393                   for (i = 1; i<np ; i++) {
394                       trackIt[i] = 1;
395                       iparticle = (TParticle *) particles->At(i);
396                       Int_t kf = iparticle->GetPdgCode();
397                       Int_t ks = iparticle->GetStatusCode();
398 // flagged particle
399
400                       if (pFlag[i] == 1) {
401                           printf("\n deselected %d", i);
402                           ipF = iparticle->GetFirstDaughter();
403                           ipL = iparticle->GetLastDaughter();   
404                           if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
405                           continue;
406                       }
407
408 // flag decay products of particles with long life-time (c tau > .3 mum)                      
409                       
410                       if (ks != 1) { 
411                           TParticlePDG *particle = pDataBase->GetParticle(kf);
412                           
413                           Double_t lifeTime = fDecayer->GetLifetime(kf);
414                           Double_t mass     = particle->Mass();
415                           Double_t width    = particle->Width();
416                           printf("\n lifetime %d %e %e %e ", i, lifeTime, mass, width);
417                           if (lifeTime > (Double_t) 1.e-15) {
418                               ipF = iparticle->GetFirstDaughter();
419                               ipL = iparticle->GetLastDaughter();       
420                               if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
421                           } else{
422                               trackIt[i]     = 0;
423                               pSelected[i]   = 1;
424                           }
425                       } // ks==1 ?
426 //
427 // children
428                       if (ChildSelected(TMath::Abs(kf)) || fForceDecay == kAll && trackIt[i])
429                       {
430                           if (fCutOnChild) {
431                               pc[0]=iparticle->Px();
432                               pc[1]=iparticle->Py();
433                               pc[2]=iparticle->Pz();
434                               Bool_t  childok = KinematicSelection(iparticle, 1);
435                               if(childok) {
436                                   pSelected[i]  = 1;
437                                   ncsel++;
438                               } else {
439                                   ncsel=-1;
440                                   break;
441                               } // child kine cuts
442                           } else {
443                               pSelected[i]  = 1;
444                               ncsel++;
445                           } // if child selection
446                       } // select muon
447                   } // decay particle loop
448               } // if decay products
449               
450               Int_t iparent;
451               if ((fCutOnChild && ncsel >0) || !fCutOnChild){
452                   ipa++;
453 //
454 // Parent
455                   SetTrack(0, -1, iPart, p, origin0, polar, 0, kPPrimary, nt, wgtp);
456                   pParent[0] = nt;
457                   KeepTrack(nt); 
458 //
459 // Decay Products
460 //                
461                   for (i = 1; i < np; i++) {
462                       if (pSelected[i]) {
463                           TParticle* iparticle = (TParticle *) particles->At(i);
464                           Int_t kf  = iparticle->GetPdgCode();
465                           Int_t ipa = iparticle->GetFirstMother()-1;
466                           
467                           och[0] = origin0[0]+iparticle->Vx()/10;
468                           och[1] = origin0[1]+iparticle->Vy()/10;
469                           och[2] = origin0[2]+iparticle->Vz()/10;
470                           pc[0]  = iparticle->Px();
471                           pc[1]  = iparticle->Py();
472                           pc[2]  = iparticle->Pz();
473                           
474                           if (ipa > -1) {
475                               iparent = pParent[ipa];
476                           } else {
477                               iparent = -1;
478                           }
479                           printf("\n SetTrack %d %d %d %d", i, kf, iparent, trackIt[i]);
480                           SetTrack(fTrackIt*trackIt[i], iparent, kf,
481                                            pc, och, polar,
482                                            0, kPDecay, nt, wgtch);
483                           pParent[i] = nt;
484                           KeepTrack(nt); 
485                       }
486                   }
487               }  // Decays by Lujet
488
489               if (pFlag)      delete[] pFlag;
490               if (pParent)    delete[] pParent;
491               if (pSelected)  delete[] pSelected;          
492               if (trackIt)    delete[] trackIt;     
493           } // kinematic selection
494           else  // nodecay option, so parent will be tracked by GEANT (pions, kaons, eta, omegas, baryons)
495           {
496             gAlice->
497                 SetTrack(fTrackIt,-1,iPart,p,origin0,polar,0,kPPrimary,nt,wgtp);
498             ipa++; 
499           }
500           break;
501     } // while
502   } // event loop
503   SetHighWaterMark(nt);
504 }
505
506 AliGenParam& AliGenParam::operator=(const  AliGenParam& rhs)
507 {
508 // Assignment operator
509     return *this;
510 }
511
512
513