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