1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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
27 <img src="picts/AliGeneratorClass.gif">
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>.
38 ///////////////////////////////////////////////////////////////////
42 #include <TClonesArray.h>
43 #include <TDatabasePDG.h>
47 #include <TParticle.h>
49 #include <TVirtualMC.h>
52 #include "AliDecayer.h"
53 #include "AliGenEventHeader.h"
54 #include "AliGenHIJINGpara.h"
58 ClassImp(AliGenHIJINGpara)
62 //_____________________________________________________________________________
63 static Double_t ptpi(const Double_t *px, const Double_t *)
66 // PT-PARAMETERIZATION CDF, PRL 61(88) 1819
67 // POWER LAW FOR PT > 500 MEV
68 // MT SCALING BELOW (T=160 MEV)
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;
79 y1 = TMath::Power(kp0 / (kp0 + kxlim), kxn);
80 xmpi2 = kxmpi * kxmpi;
81 ynorm = kb * (TMath::Exp(-sqrt(kxlim * kxlim + xmpi2) / kt ));
84 y = a * TMath::Power(kp0 / (kp0 + x), kxn);
86 y = kb* TMath::Exp(-sqrt(x * x + xmpi2) / kt);
91 //_____________________________________________________________________________
92 static Double_t ptscal(Double_t pt, Int_t np)
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};
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];
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)/
108 return fmtscal*ptpion;
111 //_____________________________________________________________________________
112 static Double_t ptka( Double_t *px, Double_t *)
115 // pt parametrisation for k
117 return ptscal(*px,2);
121 //_____________________________________________________________________________
122 static Double_t etapic( Double_t *py, Double_t *)
125 // eta parametrisation for pi
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);
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);
140 //_____________________________________________________________________________
141 static Double_t etakac( Double_t *py, Double_t *)
144 // eta parametrisation for ka
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);
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);
159 //_____________________________________________________________________________
160 AliGenHIJINGpara::AliGenHIJINGpara()
174 // Default constructor
180 //_____________________________________________________________________________
181 AliGenHIJINGpara::AliGenHIJINGpara(Int_t npart)
182 :AliGenerator(npart),
195 // Standard constructor
198 fTitle="HIJING Parametrisation Particle Generator";
203 //_____________________________________________________________________________
204 AliGenHIJINGpara::~AliGenHIJINGpara()
207 // Standard destructor
215 //_____________________________________________________________________________
216 void AliGenHIJINGpara::Init()
219 // Initialise the HIJING parametrisation
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);
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);
236 TF1 etaPic0("etaPic0",&etapic,-7,7,0);
237 TF1 etaKac0("etaKac0",&etakac,-7,7,0);
239 TF1 ptPic0("ptPic0",&ptpi,0.,15.,0);
240 TF1 ptKac0("ptKac0",&ptka,0.,15.,0);
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);
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;
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();
261 fParentWeight = (intETASel*ptFrac*phiFrac) / Float_t(fNpart);
264 fPtWgtPi = (fPtMax - fPtMin) / fPtpi->Integral(0., 20.);
265 fPtWgtKa = (fPtMax - fPtMin) / fPtka->Integral(0., 20.);
266 fParentWeight = (intETASel*phiFrac) / Float_t(fNpart);
270 AliInfo(Form("The number of particles in the selected kinematic region corresponds to %f percent of a full event",
271 100./ fParentWeight));
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));
282 if (fPi0Decays && TVirtualMC::GetMC())
283 fDecayer = TVirtualMC::GetMC()->GetDecayer();
287 fDecayer->SetForceDecay(kNeutralPion);
294 //_____________________________________________________________________________
295 void AliGenHIJINGpara::Generate()
298 // Generate one trigger
302 const Float_t kRaKpic=0.14;
303 const Float_t kBorne=1/(1+kRaKpic);
304 Float_t polar[3]= {0,0,0};
306 const Int_t kPions[3] = {kPi0, kPiPlus, kPiMinus};
307 const Int_t kKaons[4] = {kK0Long, kK0Short, kKPlus, kKMinus};
311 Float_t pt, pl, ptot, wgt;
321 for (j=0;j<3;j++) origin[j]=fOrigin[j];
324 if(fVertexSmear == kPerEvent) {
326 for (j=0; j < 3; j++) origin[j] = fVertex[j];
331 eventVertex[0] = origin[0];
332 eventVertex[1] = origin[1];
333 eventVertex[2] = origin[2];
334 Float_t eventTime = time;
336 for(i=0;i<fNpart;i++) {
339 if(random[0]<kBorne) {
340 part=kPions[Int_t (random[1]*3)];
345 part=kKaons[Int_t (random[1]*4)];
350 phi=fPhiMin+random[2]*(fPhiMax-fPhiMin);
351 theta=2*TMath::ATan(TMath::Exp(-etaf->GetRandom()));
352 if(theta<fThetaMin || theta>fThetaMax) continue;
355 pt = ptf->GetRandom();
357 pt = fPtMin + random[3] * (fPtMax - fPtMin);
361 pl=pt/TMath::Tan(theta);
362 ptot=TMath::Sqrt(pt*pt+pl*pl);
363 if(ptot<fPMin || ptot>fPMax) continue;
364 p[0]=pt*TMath::Cos(phi);
365 p[1]=pt*TMath::Sin(phi);
367 if(fVertexSmear==kPerTrack) {
370 origin[j]=fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
371 TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
374 time = fTimeOrigin + fOsigma[2]/TMath::Ccgs()*
375 TMath::Cos(2*random[0]*TMath::Pi())*
376 TMath::Sqrt(-2*TMath::Log(random[1]));
382 wgt *= (fParentWeight * ptf->Eval(pt));
386 if (part == kPi0 && fPi0Decays){
388 // Decay pi0 if requested
389 PushTrack(0,-1,part,p,origin,polar,time,kPPrimary,fNt,wgt);
391 DecayPi0(origin, p, time);
393 // printf("fNt %d", fNt);
394 PushTrack(fTrackIt,-1,part,p,origin,polar,time,kPPrimary,fNt,wgt);
401 SetHighWaterMark(fNt);
406 AliGenEventHeader* header = new AliGenEventHeader("HIJINGparam");
408 header->SetPrimaryVertex(eventVertex);
409 header->SetInteractionTime(eventTime);
410 header->SetNProduced(fNpartProd);
412 header->SetName(fName);
413 fContainer->AddHeader(header);
415 gAlice->SetGenEventHeader(header);
419 void AliGenHIJINGpara::SetPtRange(Float_t ptmin, Float_t ptmax) {
420 AliGenerator::SetPtRange(ptmin, ptmax);
423 void AliGenHIJINGpara::DecayPi0(Float_t* orig, Float_t * p, Float_t time)
427 // and put decay products on the stack
429 static TClonesArray *particles;
430 if(!particles) particles = new TClonesArray("TParticle",1000);
432 const Float_t kMass = TDatabasePDG::Instance()->GetParticle(kPi0)->Mass();
433 Float_t e = TMath::Sqrt(p[0] * p[0] + p[1] * p[1] + p[2] * p[2]+ kMass * kMass);
436 TLorentzVector pmom(p[0], p[1], p[2], e);
437 fDecayer->Decay(kPi0, &pmom);
440 // Put decay particles on the stack
442 Float_t polar[3] = {0., 0., 0.};
443 Int_t np = fDecayer->ImportParticles(particles);
444 fNpartProd += (np-1);
446 for (Int_t i = 1; i < np; i++)
448 TParticle* iParticle = (TParticle *) particles->At(i);
449 p[0] = iParticle->Px();
450 p[1] = iParticle->Py();
451 p[2] = iParticle->Pz();
452 Int_t part = iParticle->GetPdgCode();
454 PushTrack(fTrackIt, fNt, part, p, orig, polar, time, kPDecay, nt, fParentWeight);
460 void AliGenHIJINGpara::Draw( const char * /*opt*/)
463 // Draw the pT and y Distributions
465 TCanvas *c0 = new TCanvas("c0","Canvas 0",400,10,600,700);
469 fPtpi->GetHistogram()->SetXTitle("p_{T} (GeV)");
472 fPtka->GetHistogram()->SetXTitle("p_{T} (GeV)");