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