]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenHIJINGpara.cxx
track matching macros from Alberto
[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
62 //_____________________________________________________________________________
63 static Double_t ptpi(Double_t *px, Double_t *)
64 {
65   //
66   //     PT-PARAMETERIZATION CDF, PRL 61(88) 1819
67   //     POWER LAW FOR PT > 500 MEV
68   //     MT SCALING BELOW (T=160 MEV)
69   //
70   const Double_t kp0    = 1.3;
71   const Double_t kxn    = 8.28;
72   const Double_t kxlim  = 0.5;
73   const Double_t kt     = 0.160;
74   const Double_t kxmpi  = 0.139;
75   const Double_t kb     = 1.;
76   Double_t y, y1, xmpi2, ynorm, a;
77   Double_t x = *px;
78   //
79   y1 = TMath::Power(kp0 / (kp0 + kxlim), kxn);
80   xmpi2 = kxmpi * kxmpi;
81   ynorm = kb * (TMath::Exp(-sqrt(kxlim * kxlim + xmpi2) / kt ));
82   a = ynorm / y1;
83   if (x > kxlim)
84     y = a * TMath::Power(kp0 / (kp0 + x), kxn);
85   else
86     y = kb* TMath::Exp(-sqrt(x * x + xmpi2) / kt);
87   
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         fNt(-1),
163         fNpartProd(0),
164         fPi0Decays(kFALSE),
165         fPtWgtPi(0.),
166         fPtWgtKa(0.),
167         fPtpi(0),
168         fPtka(0),
169         fETApic(0),
170         fETAkac(0),
171         fDecayer(0)
172 {
173     //
174     // Default constructor
175     //
176     SetCutVertexZ();
177     SetPtRange();
178 }
179
180 //_____________________________________________________________________________
181 AliGenHIJINGpara::AliGenHIJINGpara(Int_t npart)
182     :AliGenerator(npart),
183         fNt(-1),
184         fNpartProd(npart),
185         fPi0Decays(kFALSE),
186         fPtWgtPi(0.),
187         fPtWgtKa(0.),
188         fPtpi(0),
189         fPtka(0),
190         fETApic(0),
191         fETAkac(0),
192         fDecayer(0)
193 {
194   // 
195   // Standard constructor
196   //
197     fName="HIJINGpara";
198     fTitle="HIJING Parametrisation Particle Generator";
199     SetCutVertexZ();
200     SetPtRange();
201 }
202
203 AliGenHIJINGpara::AliGenHIJINGpara(const AliGenHIJINGpara & para):
204     AliGenerator(para),
205     fNt(-1),
206     fNpartProd(0),
207     fPi0Decays(kFALSE),
208     fPtWgtPi(0.),
209     fPtWgtKa(0.),
210     fPtpi(0),
211     fPtka(0),
212     fETApic(0),
213     fETAkac(0),
214     fDecayer(0)
215 {
216 // Copy constructor
217     para.Copy(*this);
218 }
219
220 //_____________________________________________________________________________
221 AliGenHIJINGpara::~AliGenHIJINGpara()
222 {
223   //
224   // Standard destructor
225   //
226     delete fPtpi;
227     delete fPtka;
228     delete fETApic;
229     delete fETAkac;
230 }
231
232 //_____________________________________________________________________________
233 void AliGenHIJINGpara::Init()
234 {
235   //
236   // Initialise the HIJING parametrisation
237   //
238     Float_t etaMin =-TMath::Log(TMath::Tan(
239         TMath::Min((Double_t)fThetaMax/2,TMath::Pi()/2-1.e-10)));
240     Float_t etaMax = -TMath::Log(TMath::Tan(
241         TMath::Max((Double_t)fThetaMin/2,1.e-10)));
242     fPtpi   = new TF1("ptpi",&ptpi,0,20,0);
243     gROOT->GetListOfFunctions()->Remove(fPtpi);
244     fPtka   = new TF1("ptka",&ptka,0,20,0);
245     gROOT->GetListOfFunctions()->Remove(fPtka);
246     fPtpi->SetNpx(1000);
247     fPtka->SetNpx(1000);
248     fETApic = new TF1("etapic",&etapic,etaMin,etaMax,0);
249     gROOT->GetListOfFunctions()->Remove(fETApic);
250     fETAkac = new TF1("etakac",&etakac,etaMin,etaMax,0);
251     gROOT->GetListOfFunctions()->Remove(fETAkac);
252
253     TF1 etaPic0("etaPic0",&etapic,-7,7,0);
254     TF1 etaKac0("etaKac0",&etakac,-7,7,0);
255
256     TF1 ptPic0("ptPic0",&ptpi,0.,15.,0);
257     TF1 ptKac0("ptKac0",&ptka,0.,15.,0);
258
259     Float_t intETApi  = etaPic0.Integral(-0.5, 0.5);
260     Float_t intETAka  = etaKac0.Integral(-0.5, 0.5);
261     Float_t scalePi   = 7316/(intETApi/1.5);
262     Float_t scaleKa   =  684/(intETAka/2.0);
263
264 //  Fraction of events corresponding to the selected pt-range    
265     Float_t intPt    = (0.877*ptPic0.Integral(0, 15)+
266                         0.123*ptKac0.Integral(0, 15));
267     Float_t intPtSel = (0.877*ptPic0.Integral(fPtMin, fPtMax)+
268                         0.123*ptKac0.Integral(fPtMin, fPtMax));
269     Float_t ptFrac   = intPtSel/intPt;
270
271 //  Fraction of events corresponding to the selected eta-range    
272     Float_t intETASel  = (scalePi*etaPic0.Integral(etaMin, etaMax)+
273                           scaleKa*etaKac0.Integral(etaMin, etaMax));
274 //  Fraction of events corresponding to the selected phi-range    
275     Float_t phiFrac    = (fPhiMax-fPhiMin)/2/TMath::Pi();
276
277     
278     fParentWeight = (intETASel*ptFrac*phiFrac) / Float_t(fNpart);
279     
280     if (fAnalog != 0) {
281         fPtWgtPi = (fPtMax - fPtMin) / fPtpi->Integral(0., 20.);
282         fPtWgtKa = (fPtMax - fPtMin) / fPtka->Integral(0., 20.);
283         fParentWeight = (intETASel*phiFrac) / Float_t(fNpart);
284     }
285     
286     
287     AliInfo(Form("The number of particles in the selected kinematic region corresponds to %f percent of a full event", 
288                  100./ fParentWeight));
289
290 // Issue warning message if etaMin or etaMax are outside the alowed range 
291 // of the parametrization
292     if (etaMin < -8.001 || etaMax > 8.001) {
293         AliWarning("\nYOU ARE USING THE PARAMETERISATION OUTSIDE ");    
294         AliWarning("THE ALLOWED PSEUDORAPIDITY RANGE (-8. - 8.)");          
295         AliWarning(Form("YOUR LIMITS: %f %f \n ", etaMin, etaMax));
296     }
297 //
298 //
299     if (fPi0Decays && gMC)
300         fDecayer = gMC->GetDecayer();
301 }
302
303
304 //_____________________________________________________________________________
305 void AliGenHIJINGpara::Generate()
306 {
307   //
308   // Generate one trigger
309   //
310
311   
312     const Float_t kRaKpic=0.14;
313     const Float_t kBorne=1/(1+kRaKpic);
314     Float_t polar[3]= {0,0,0};
315     //
316     const Int_t kPions[3] = {kPi0, kPiPlus, kPiMinus};
317     const Int_t kKaons[4] = {kK0Long, kK0Short, kKPlus, kKMinus};
318     //
319     Float_t origin[3];
320     Float_t pt, pl, ptot, wgt;
321     Float_t phi, theta;
322     Float_t p[3];
323     Int_t i, part, j;
324     //
325     TF1 *ptf;
326     TF1 *etaf;
327     //
328     Float_t random[6];
329     //
330     for (j=0;j<3;j++) origin[j]=fOrigin[j];
331
332     if(fVertexSmear == kPerEvent) {
333         Vertex();
334         for (j=0; j < 3; j++) origin[j] = fVertex[j];
335     } // if kPerEvent
336     TArrayF eventVertex;
337     eventVertex.Set(3);
338     eventVertex[0] = origin[0];
339     eventVertex[1] = origin[1];
340     eventVertex[2] = origin[2];
341      
342     for(i=0;i<fNpart;i++) {
343         while(1) {
344             Rndm(random,4);
345             if(random[0]<kBorne) {
346                 part=kPions[Int_t (random[1]*3)];
347                 ptf=fPtpi;
348                 etaf=fETApic;
349                 wgt = fPtWgtPi;
350             } else {
351                 part=kKaons[Int_t (random[1]*4)];
352                 ptf=fPtka;
353                 etaf=fETAkac;
354                 wgt = fPtWgtKa;
355             }
356             phi=fPhiMin+random[2]*(fPhiMax-fPhiMin);
357             theta=2*TMath::ATan(TMath::Exp(-etaf->GetRandom()));
358             if(theta<fThetaMin || theta>fThetaMax) continue;
359          
360             if (fAnalog == 0) { 
361                 pt = ptf->GetRandom();
362             } else {
363                 pt = fPtMin + random[3] * (fPtMax - fPtMin);
364             }
365             
366             
367             pl=pt/TMath::Tan(theta);
368             ptot=TMath::Sqrt(pt*pt+pl*pl);
369             if(ptot<fPMin || ptot>fPMax) continue;
370             p[0]=pt*TMath::Cos(phi);
371             p[1]=pt*TMath::Sin(phi);
372             p[2]=pl;
373             if(fVertexSmear==kPerTrack) {
374                 Rndm(random,6);
375                 for (j=0;j<3;j++) {
376                     origin[j]=fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
377                         TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
378                 }
379             }
380
381             if (fAnalog == 0) { 
382                 wgt = fParentWeight;
383             } else {
384                 wgt *= (fParentWeight * ptf->Eval(pt));
385             }
386             
387             
388             if (part == kPi0 && fPi0Decays){
389 //
390 //          Decay pi0 if requested
391                 PushTrack(0,-1,part,p,origin,polar,0,kPPrimary,fNt,wgt);
392                 KeepTrack(fNt);
393                 DecayPi0(origin, p);
394             } else {
395       // printf("fNt %d", fNt);
396                 PushTrack(fTrackIt,-1,part,p,origin,polar,0,kPPrimary,fNt,wgt);
397
398                 KeepTrack(fNt);
399             }
400
401             break;
402         }
403         SetHighWaterMark(fNt);
404     }
405 //
406
407 // Header
408     AliGenEventHeader* header = new AliGenEventHeader("HIJINGparam");
409 // Event Vertex
410     header->SetPrimaryVertex(eventVertex);
411     header->SetNProduced(fNpartProd);
412     gAlice->SetGenEventHeader(header); 
413 }
414
415 void AliGenHIJINGpara::SetPtRange(Float_t ptmin, Float_t ptmax) {
416     AliGenerator::SetPtRange(ptmin, ptmax);
417 }
418
419 void AliGenHIJINGpara::DecayPi0(Float_t* orig, Float_t * p) 
420 {
421 //
422 //    Decay the pi0
423 //    and put decay products on the stack
424 //
425     static TClonesArray *particles;
426     if(!particles) particles = new TClonesArray("TParticle",1000);
427 //    
428     const Float_t kMass = TDatabasePDG::Instance()->GetParticle(kPi0)->Mass();
429     Float_t       e     = TMath::Sqrt(p[0] * p[0] + p[1] * p[1] + p[2] * p[2]+ kMass * kMass);
430 //
431 //  Decay the pi0    
432     TLorentzVector pmom(p[0], p[1], p[2], e);
433     fDecayer->Decay(kPi0, &pmom);
434     
435 //
436 // Put decay particles on the stack
437 //
438     Float_t polar[3] = {0., 0., 0.};
439     Int_t np = fDecayer->ImportParticles(particles);
440     fNpartProd += (np-1);
441     Int_t nt;    
442     for (Int_t i = 1; i < np; i++)
443     {
444         TParticle* iParticle =  (TParticle *) particles->At(i);
445         p[0] = iParticle->Px();
446         p[1] = iParticle->Py();
447         p[2] = iParticle->Pz();
448         Int_t part = iParticle->GetPdgCode();
449
450         PushTrack(fTrackIt, fNt, part, p, orig, polar, 0, kPDecay, nt, fParentWeight);
451         KeepTrack(nt);
452     }
453     fNt = nt;
454 }
455
456 void AliGenHIJINGpara::Copy(TObject &) const
457 {
458   Fatal("Copy","Not implemented!\n");
459 }
460
461
462 void AliGenHIJINGpara::Draw( const char * /*opt*/)
463 {
464     //
465     // Draw the pT and y Distributions
466     //
467      TCanvas *c0 = new TCanvas("c0","Canvas 0",400,10,600,700);
468      c0->Divide(2,1);
469      c0->cd(1);
470      fPtpi->Draw();
471      fPtpi->GetHistogram()->SetXTitle("p_{T} (GeV)");     
472      c0->cd(2);
473      fPtka->Draw();
474      fPtka->GetHistogram()->SetXTitle("p_{T} (GeV)");     
475
476 }