]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenHIJINGpara.cxx
Change to avoid name conflict with Root v3-10-01a
[u/mrichter/AliRoot.git] / EVGEN / AliGenHIJINGpara.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 // Parameterisation of pi and K, eta and pt distributions
19 // used for the ALICE TDRs.
20 // eta: according to HIJING (shadowing + quenching)
21 // pT : according to CDF measurement at 1.8 TeV
22 // Author: andreas.morsch@cern.ch
23
24
25 //Begin_Html
26 /*
27 <img src="picts/AliGeneratorClass.gif">
28 </pre>
29 <br clear=left>
30 <font size=+2 color=red>
31 <p>The responsible person for this module is
32 <a href="mailto:andreas.morsch@cern.ch">Andreas Morsch</a>.
33 </font>
34 <pre>
35 */
36 //End_Html
37 //                                                               //
38 ///////////////////////////////////////////////////////////////////
39
40 #include <TArrayF.h>
41 #include <TClonesArray.h>
42 #include <TDatabasePDG.h>
43 #include <TF1.h>
44 #include <TParticle.h>
45 #include <TPDGCode.h>
46 #include <TVirtualMC.h>
47
48 #include "AliConst.h"
49 #include "AliDecayer.h"
50 #include "AliGenEventHeader.h"
51 #include "AliGenHIJINGpara.h"
52 #include "AliRun.h"
53
54 ClassImp(AliGenHIJINGpara)
55
56
57 AliGenHIJINGpara::AliGenHIJINGpara(const AliGenHIJINGpara & para):
58     AliGenerator(para)
59 {
60 // Copy constructor
61     para.Copy(*this);
62 }
63
64 //_____________________________________________________________________________
65 static Double_t ptpi(Double_t *px, Double_t *)
66 {
67   //
68   //     PT-PARAMETERIZATION CDF, PRL 61(88) 1819
69   //     POWER LAW FOR PT > 500 MEV
70   //     MT SCALING BELOW (T=160 MEV)
71   //
72   const Double_t kp0 = 1.3;
73   const Double_t kxn = 8.28;
74   const Double_t kxlim=0.5;
75   const Double_t kt=0.160;
76   const Double_t kxmpi=0.139;
77   const Double_t kb=1.;
78   Double_t y, y1, xmpi2, ynorm, a;
79   Double_t x=*px;
80   //
81   y1=TMath::Power(kp0/(kp0+kxlim),kxn);
82   xmpi2=kxmpi*kxmpi;
83   ynorm=kb*(TMath::Exp(-sqrt(kxlim*kxlim+xmpi2)/kt));
84   a=ynorm/y1;
85   if (x > kxlim)
86     y=a*TMath::Power(kp0/(kp0+x),kxn);
87   else
88     y=kb*TMath::Exp(-sqrt(x*x+xmpi2)/kt);
89   return y*x;
90 }
91
92 //_____________________________________________________________________________
93 static Double_t ptscal(Double_t pt, Int_t np)
94 {
95     //    SCALING EN MASSE PAR RAPPORT A PTPI
96     //     MASS PI,K,ETA,RHO,OMEGA,ETA',PHI
97     const Double_t khm[10] = {.13957,.493,.5488,.769,.7826,.958,1.02,0,0,0};
98     //     VALUE MESON/PI AT 5 GEV
99     const Double_t kfmax[10]={1.,0.3,0.55,1.0,1.0,1.0,1.0,0,0,0};
100     np--;
101     Double_t f5=TMath::Power(((
102         sqrt(100.018215)+2.)/(sqrt(100.+khm[np]*khm[np])+2.0)),12.3);
103     Double_t fmax2=f5/kfmax[np];
104     // PIONS
105     Double_t ptpion=100.*ptpi(&pt, (Double_t*) 0);
106     Double_t fmtscal=TMath::Power(((
107         sqrt(pt*pt+0.018215)+2.)/ (sqrt(pt*pt+khm[np]*khm[np])+2.0)),12.3)/ 
108         fmax2;
109     return fmtscal*ptpion;
110 }
111
112 //_____________________________________________________________________________
113 static Double_t ptka( Double_t *px, Double_t *)
114 {
115     //
116     // pt parametrisation for k
117     //
118     return ptscal(*px,2);
119 }
120
121
122 //_____________________________________________________________________________
123 static Double_t etapic( Double_t *py, Double_t *)
124 {
125   //
126   // eta parametrisation for pi
127   //
128     const Double_t ka1    = 4913.;
129     const Double_t ka2    = 1819.;
130     const Double_t keta1  = 0.22;
131     const Double_t keta2  = 3.66;
132     const Double_t kdeta1 = 1.47;
133     const Double_t kdeta2 = 1.51;
134     Double_t y=TMath::Abs(*py);
135     //
136     Double_t ex1 = (y-keta1)*(y-keta1)/(2*kdeta1*kdeta1);
137     Double_t ex2 = (y-keta2)*(y-keta2)/(2*kdeta2*kdeta2);
138     return ka1*TMath::Exp(-ex1)+ka2*TMath::Exp(-ex2);
139 }
140
141 //_____________________________________________________________________________
142 static Double_t etakac( Double_t *py, Double_t *)
143 {
144     //
145     // eta parametrisation for ka
146     //
147     const Double_t ka1    = 497.6;
148     const Double_t ka2    = 215.6;
149     const Double_t keta1  = 0.79;
150     const Double_t keta2  = 4.09;
151     const Double_t kdeta1 = 1.54;
152     const Double_t kdeta2 = 1.40;
153     Double_t y=TMath::Abs(*py);
154     //
155     Double_t ex1 = (y-keta1)*(y-keta1)/(2*kdeta1*kdeta1);
156     Double_t ex2 = (y-keta2)*(y-keta2)/(2*kdeta2*kdeta2);
157     return ka1*TMath::Exp(-ex1)+ka2*TMath::Exp(-ex2);
158 }
159
160 //_____________________________________________________________________________
161 AliGenHIJINGpara::AliGenHIJINGpara()
162   :AliGenerator()
163 {
164     //
165     // Default constructor
166     //
167     fPtpi    =  0;
168     fPtka    =  0;
169     fETApic  =  0;
170     fETAkac  =  0;
171     fDecayer =  0;
172     fNt      = -1;
173     SetCutVertexZ();
174     SetPtRange();
175     SetPi0Decays();
176 }
177
178 //_____________________________________________________________________________
179 AliGenHIJINGpara::AliGenHIJINGpara(Int_t npart)
180   :AliGenerator(npart)
181 {
182   // 
183   // Standard constructor
184   //
185     fName="HIJINGpara";
186     fTitle="HIJING Parametrisation Particle Generator";
187     fPtpi    =  0;
188     fPtka    =  0;
189     fETApic  =  0;
190     fETAkac  =  0;
191     fDecayer =  0;
192     fNt      = -1;
193     SetCutVertexZ();
194     SetPtRange();
195     SetPi0Decays();
196 }
197
198 //_____________________________________________________________________________
199 AliGenHIJINGpara::~AliGenHIJINGpara()
200 {
201   //
202   // Standard destructor
203   //
204     delete fPtpi;
205     delete fPtka;
206     delete fETApic;
207     delete fETAkac;
208 }
209
210 //_____________________________________________________________________________
211 void AliGenHIJINGpara::Init()
212 {
213   //
214   // Initialise the HIJING parametrisation
215   //
216     Float_t etaMin =-TMath::Log(TMath::Tan(
217         TMath::Min((Double_t)fThetaMax/2,TMath::Pi()/2-1.e-10)));
218     Float_t etaMax = -TMath::Log(TMath::Tan(
219         TMath::Max((Double_t)fThetaMin/2,1.e-10)));
220     fPtpi   = new TF1("ptpi",&ptpi,0,20,0);
221     fPtka   = new TF1("ptka",&ptka,0,20,0);
222     fETApic = new TF1("etapic",&etapic,etaMin,etaMax,0);
223     fETAkac = new TF1("etakac",&etakac,etaMin,etaMax,0);
224
225     TF1 etaPic0("etaPic0",&etapic,-7,7,0);
226     TF1 etaKac0("etaKac0",&etakac,-7,7,0);
227
228     TF1 ptPic0("ptPic0",&ptpi,0.,15.,0);
229     TF1 ptKac0("ptKac0",&ptka,0.,15.,0);
230
231     Float_t intETApi  = etaPic0.Integral(-0.5, 0.5);
232     Float_t intETAka  = etaKac0.Integral(-0.5, 0.5);
233     Float_t scalePi   = 7316/(intETApi/1.5);
234     Float_t scaleKa   =  684/(intETAka/2.0);
235
236 //  Fraction of events corresponding to the selected pt-range    
237     Float_t intPt    = (0.877*ptPic0.Integral(0, 15)+
238                         0.123*ptKac0.Integral(0, 15));
239     Float_t intPtSel = (0.877*ptPic0.Integral(fPtMin, fPtMax)+
240                         0.123*ptKac0.Integral(fPtMin, fPtMax));
241     Float_t ptFrac   = intPtSel/intPt;
242
243 //  Fraction of events corresponding to the selected eta-range    
244     Float_t intETASel  = (scalePi*etaPic0.Integral(etaMin, etaMax)+
245                           scaleKa*etaKac0.Integral(etaMin, etaMax));
246 //  Fraction of events corresponding to the selected phi-range    
247     Float_t phiFrac    = (fPhiMax-fPhiMin)/2/TMath::Pi();
248
249     
250     fParentWeight = Float_t(fNpart)/(intETASel*ptFrac*phiFrac);
251     
252     if (fAnalog != 0) {
253         fPtWgtPi = (fPtMax - fPtMin) / fPtpi->Integral(0., 20.);
254         fPtWgtKa = (fPtMax - fPtMin) / fPtka->Integral(0., 20.);
255         fParentWeight = Float_t(fNpart)/(intETASel*phiFrac);
256     }
257     
258     
259     printf("%s: The number of particles in the selected kinematic region corresponds to %f percent of a full event\n ", 
260            ClassName(),100.*fParentWeight);
261
262 // Issue warning message if etaMin or etaMax are outside the alowed range 
263 // of the parametrization
264     if (etaMin < -8.001 || etaMax > 8.001) {
265         printf("\n \n WARNING FROM AliGenHIJINGPara !");
266         printf("\n YOU ARE USING THE PARAMETERISATION OUTSIDE ");       
267         printf("\n THE ALLOWED PSEUDORAPIDITY RANGE (-8. - 8.)");           
268         printf("\n YOUR LIMITS: %f %f \n \n ", etaMin, etaMax);
269     }
270 //
271 //
272     if (fPi0Decays && gMC)
273         fDecayer = gMC->GetDecayer();
274 }
275
276
277 //_____________________________________________________________________________
278 void AliGenHIJINGpara::Generate()
279 {
280   //
281   // Generate one trigger
282   //
283
284   
285     const Float_t kRaKpic=0.14;
286     const Float_t kBorne=1/(1+kRaKpic);
287     Float_t polar[3]= {0,0,0};
288     //
289     const Int_t kPions[3] = {kPi0, kPiPlus, kPiMinus};
290     const Int_t kKaons[4] = {kK0Long, kK0Short, kKPlus, kKMinus};
291     //
292     Float_t origin[3];
293     Float_t pt, pl, ptot, wgt;
294     Float_t phi, theta;
295     Float_t p[3];
296     Int_t i, part, j;
297     //
298     TF1 *ptf;
299     TF1 *etaf;
300     //
301     Float_t random[6];
302     //
303     for (j=0;j<3;j++) origin[j]=fOrigin[j];
304
305     if(fVertexSmear == kPerEvent) {
306         Vertex();
307         for (j=0; j < 3; j++) origin[j] = fVertex[j];
308     } // if kPerEvent
309     TArrayF eventVertex;
310     eventVertex.Set(3);
311     eventVertex[0] = origin[0];
312     eventVertex[1] = origin[1];
313     eventVertex[2] = origin[2];
314      
315     for(i=0;i<fNpart;i++) {
316         while(1) {
317             Rndm(random,4);
318             if(random[0]<kBorne) {
319                 part=kPions[Int_t (random[1]*3)];
320                 ptf=fPtpi;
321                 etaf=fETApic;
322                 wgt = fPtWgtPi;
323             } else {
324                 part=kKaons[Int_t (random[1]*4)];
325                 ptf=fPtka;
326                 etaf=fETAkac;
327                 wgt = fPtWgtKa;
328             }
329             phi=fPhiMin+random[2]*(fPhiMax-fPhiMin);
330             theta=2*TMath::ATan(TMath::Exp(-etaf->GetRandom()));
331             if(theta<fThetaMin || theta>fThetaMax) continue;
332          
333             if (fAnalog == 0) { 
334                 pt = ptf->GetRandom();
335             } else {
336                 pt = fPtMin + random[3] * (fPtMax - fPtMin);
337             }
338             
339             
340             pl=pt/TMath::Tan(theta);
341             ptot=TMath::Sqrt(pt*pt+pl*pl);
342             if(ptot<fPMin || ptot>fPMax) continue;
343             p[0]=pt*TMath::Cos(phi);
344             p[1]=pt*TMath::Sin(phi);
345             p[2]=pl;
346             if(fVertexSmear==kPerTrack) {
347                 Rndm(random,6);
348                 for (j=0;j<3;j++) {
349                     origin[j]=fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
350                         TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
351                 }
352             }
353             
354             if (fAnalog == 0) { 
355                 wgt = fParentWeight;
356             } else {
357                 wgt *= (fParentWeight * ptf->Eval(pt));
358             }
359             
360             
361             if (part == kPi0 && fPi0Decays){
362 //
363 //          Decay pi0 if requested
364                 PushTrack(0,-1,part,p,origin,polar,0,kPPrimary,fNt,fParentWeight);
365                 KeepTrack(fNt);
366                 DecayPi0(origin, p);
367             } else {
368                 PushTrack(fTrackIt,-1,part,p,origin,polar,0,kPPrimary,fNt,fParentWeight);
369                 KeepTrack(fNt);
370             }
371
372             break;
373         }
374         SetHighWaterMark(fNt);
375     }
376 //
377
378 // Header
379     AliGenEventHeader* header = new AliGenEventHeader("HIJINGparam");
380 // Event Vertex
381     header->SetPrimaryVertex(eventVertex);
382     gAlice->SetGenEventHeader(header); 
383 }
384
385 void AliGenHIJINGpara::SetPtRange(Float_t ptmin, Float_t ptmax) {
386     AliGenerator::SetPtRange(ptmin, ptmax);
387 }
388
389 void AliGenHIJINGpara::DecayPi0(Float_t* orig, Float_t * p) 
390 {
391 //
392 //    Decay the pi0
393 //    and put decay products on the stack
394 //
395     static TClonesArray *particles;
396     if(!particles) particles = new TClonesArray("TParticle",1000);
397 //    
398     const Float_t kMass = TDatabasePDG::Instance()->GetParticle(kPi0)->Mass();
399     Float_t       e     = TMath::Sqrt(p[0] * p[0] + p[1] * p[1] + p[2] * p[2]+ kMass * kMass);
400 //
401 //  Decay the pi0    
402     TLorentzVector pmom(p[0], p[1], p[2], e);
403     fDecayer->Decay(kPi0, &pmom);
404     
405 //
406 // Put decay particles on the stack
407 //
408     Float_t polar[3] = {0., 0., 0.};
409     Int_t np = fDecayer->ImportParticles(particles);
410     Int_t nt;    
411     for (Int_t i = 1; i < np; i++)
412     {
413         TParticle* iParticle =  (TParticle *) particles->At(i);
414         p[0] = iParticle->Px();
415         p[1] = iParticle->Py();
416         p[2] = iParticle->Pz();
417         Int_t part = iParticle->GetPdgCode();
418
419         PushTrack(fTrackIt, fNt, part, p, orig, polar, 0, kPDecay, nt, fParentWeight);
420         KeepTrack(nt);
421     }
422     fNt = nt;
423 }
424
425 void AliGenHIJINGpara::Copy(AliGenHIJINGpara &) const
426 {
427   Fatal("Copy","Not implemented!\n");
428 }