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