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