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)
61 AliGenHIJINGpara::AliGenHIJINGpara(const AliGenHIJINGpara & para):
68 //_____________________________________________________________________________
69 static Double_t ptpi(Double_t *px, Double_t *)
72 // PT-PARAMETERIZATION CDF, PRL 61(88) 1819
73 // POWER LAW FOR PT > 500 MEV
74 // MT SCALING BELOW (T=160 MEV)
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;
85 y1 = TMath::Power(kp0 / (kp0 + kxlim), kxn);
86 xmpi2 = kxmpi * kxmpi;
87 ynorm = kb * (TMath::Exp(-sqrt(kxlim * kxlim + xmpi2) / kt ));
90 y = a * TMath::Power(kp0 / (kp0 + x), kxn);
92 y = kb* TMath::Exp(-sqrt(x * x + xmpi2) / kt);
97 //_____________________________________________________________________________
98 static Double_t ptscal(Double_t pt, Int_t np)
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};
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];
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)/
114 return fmtscal*ptpion;
117 //_____________________________________________________________________________
118 static Double_t ptka( Double_t *px, Double_t *)
121 // pt parametrisation for k
123 return ptscal(*px,2);
127 //_____________________________________________________________________________
128 static Double_t etapic( Double_t *py, Double_t *)
131 // eta parametrisation for pi
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);
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);
146 //_____________________________________________________________________________
147 static Double_t etakac( Double_t *py, Double_t *)
150 // eta parametrisation for ka
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);
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);
165 //_____________________________________________________________________________
166 AliGenHIJINGpara::AliGenHIJINGpara()
170 // Default constructor
184 //_____________________________________________________________________________
185 AliGenHIJINGpara::AliGenHIJINGpara(Int_t npart)
189 // Standard constructor
192 fTitle="HIJING Parametrisation Particle Generator";
205 //_____________________________________________________________________________
206 AliGenHIJINGpara::~AliGenHIJINGpara()
209 // Standard destructor
217 //_____________________________________________________________________________
218 void AliGenHIJINGpara::Init()
221 // Initialise the HIJING parametrisation
223 Float_t etaMin =-TMath::Log(TMath::Tan(
224 TMath::Min((Double_t)fThetaMax/2,TMath::Pi()/2-1.e-10)));
225 Float_t etaMax = -TMath::Log(TMath::Tan(
226 TMath::Max((Double_t)fThetaMin/2,1.e-10)));
227 fPtpi = new TF1("ptpi",&ptpi,0,20,0);
228 gROOT->GetListOfFunctions()->Remove(fPtpi);
229 fPtka = new TF1("ptka",&ptka,0,20,0);
230 gROOT->GetListOfFunctions()->Remove(fPtka);
233 fETApic = new TF1("etapic",&etapic,etaMin,etaMax,0);
234 gROOT->GetListOfFunctions()->Remove(fETApic);
235 fETAkac = new TF1("etakac",&etakac,etaMin,etaMax,0);
236 gROOT->GetListOfFunctions()->Remove(fETAkac);
238 TF1 etaPic0("etaPic0",&etapic,-7,7,0);
239 TF1 etaKac0("etaKac0",&etakac,-7,7,0);
241 TF1 ptPic0("ptPic0",&ptpi,0.,15.,0);
242 TF1 ptKac0("ptKac0",&ptka,0.,15.,0);
244 Float_t intETApi = etaPic0.Integral(-0.5, 0.5);
245 Float_t intETAka = etaKac0.Integral(-0.5, 0.5);
246 Float_t scalePi = 7316/(intETApi/1.5);
247 Float_t scaleKa = 684/(intETAka/2.0);
249 // Fraction of events corresponding to the selected pt-range
250 Float_t intPt = (0.877*ptPic0.Integral(0, 15)+
251 0.123*ptKac0.Integral(0, 15));
252 Float_t intPtSel = (0.877*ptPic0.Integral(fPtMin, fPtMax)+
253 0.123*ptKac0.Integral(fPtMin, fPtMax));
254 Float_t ptFrac = intPtSel/intPt;
256 // Fraction of events corresponding to the selected eta-range
257 Float_t intETASel = (scalePi*etaPic0.Integral(etaMin, etaMax)+
258 scaleKa*etaKac0.Integral(etaMin, etaMax));
259 // Fraction of events corresponding to the selected phi-range
260 Float_t phiFrac = (fPhiMax-fPhiMin)/2/TMath::Pi();
263 fParentWeight = (intETASel*ptFrac*phiFrac) / Float_t(fNpart);
266 fPtWgtPi = (fPtMax - fPtMin) / fPtpi->Integral(0., 20.);
267 fPtWgtKa = (fPtMax - fPtMin) / fPtka->Integral(0., 20.);
268 fParentWeight = (intETASel*phiFrac) / Float_t(fNpart);
272 AliInfo(Form("The number of particles in the selected kinematic region corresponds to %f percent of a full event",
273 100./ fParentWeight));
275 // Issue warning message if etaMin or etaMax are outside the alowed range
276 // of the parametrization
277 if (etaMin < -8.001 || etaMax > 8.001) {
278 AliWarning("\nYOU ARE USING THE PARAMETERISATION OUTSIDE ");
279 AliWarning("THE ALLOWED PSEUDORAPIDITY RANGE (-8. - 8.)");
280 AliWarning(Form("YOUR LIMITS: %f %f \n ", etaMin, etaMax));
284 if (fPi0Decays && gMC)
285 fDecayer = gMC->GetDecayer();
289 //_____________________________________________________________________________
290 void AliGenHIJINGpara::Generate()
293 // Generate one trigger
297 const Float_t kRaKpic=0.14;
298 const Float_t kBorne=1/(1+kRaKpic);
299 Float_t polar[3]= {0,0,0};
301 const Int_t kPions[3] = {kPi0, kPiPlus, kPiMinus};
302 const Int_t kKaons[4] = {kK0Long, kK0Short, kKPlus, kKMinus};
305 Float_t pt, pl, ptot, wgt;
315 for (j=0;j<3;j++) origin[j]=fOrigin[j];
317 if(fVertexSmear == kPerEvent) {
319 for (j=0; j < 3; j++) origin[j] = fVertex[j];
323 eventVertex[0] = origin[0];
324 eventVertex[1] = origin[1];
325 eventVertex[2] = origin[2];
327 for(i=0;i<fNpart;i++) {
330 if(random[0]<kBorne) {
331 part=kPions[Int_t (random[1]*3)];
336 part=kKaons[Int_t (random[1]*4)];
341 phi=fPhiMin+random[2]*(fPhiMax-fPhiMin);
342 theta=2*TMath::ATan(TMath::Exp(-etaf->GetRandom()));
343 if(theta<fThetaMin || theta>fThetaMax) continue;
346 pt = ptf->GetRandom();
348 pt = fPtMin + random[3] * (fPtMax - fPtMin);
352 pl=pt/TMath::Tan(theta);
353 ptot=TMath::Sqrt(pt*pt+pl*pl);
354 if(ptot<fPMin || ptot>fPMax) continue;
355 p[0]=pt*TMath::Cos(phi);
356 p[1]=pt*TMath::Sin(phi);
358 if(fVertexSmear==kPerTrack) {
361 origin[j]=fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
362 TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
369 wgt *= (fParentWeight * ptf->Eval(pt));
373 if (part == kPi0 && fPi0Decays){
375 // Decay pi0 if requested
376 PushTrack(0,-1,part,p,origin,polar,0,kPPrimary,fNt,wgt);
380 // printf("fNt %d", fNt);
381 PushTrack(fTrackIt,-1,part,p,origin,polar,0,kPPrimary,fNt,wgt);
388 SetHighWaterMark(fNt);
393 AliGenEventHeader* header = new AliGenEventHeader("HIJINGparam");
395 header->SetPrimaryVertex(eventVertex);
396 header->SetNProduced(fNpartProd);
397 gAlice->SetGenEventHeader(header);
400 void AliGenHIJINGpara::SetPtRange(Float_t ptmin, Float_t ptmax) {
401 AliGenerator::SetPtRange(ptmin, ptmax);
404 void AliGenHIJINGpara::DecayPi0(Float_t* orig, Float_t * p)
408 // and put decay products on the stack
410 static TClonesArray *particles;
411 if(!particles) particles = new TClonesArray("TParticle",1000);
413 const Float_t kMass = TDatabasePDG::Instance()->GetParticle(kPi0)->Mass();
414 Float_t e = TMath::Sqrt(p[0] * p[0] + p[1] * p[1] + p[2] * p[2]+ kMass * kMass);
417 TLorentzVector pmom(p[0], p[1], p[2], e);
418 fDecayer->Decay(kPi0, &pmom);
421 // Put decay particles on the stack
423 Float_t polar[3] = {0., 0., 0.};
424 Int_t np = fDecayer->ImportParticles(particles);
425 fNpartProd += (np-1);
427 for (Int_t i = 1; i < np; i++)
429 TParticle* iParticle = (TParticle *) particles->At(i);
430 p[0] = iParticle->Px();
431 p[1] = iParticle->Py();
432 p[2] = iParticle->Pz();
433 Int_t part = iParticle->GetPdgCode();
435 PushTrack(fTrackIt, fNt, part, p, orig, polar, 0, kPDecay, nt, fParentWeight);
441 void AliGenHIJINGpara::Copy(TObject &) const
443 Fatal("Copy","Not implemented!\n");
447 void AliGenHIJINGpara::Draw( const char * /*opt*/)
450 // Draw the pT and y Distributions
452 TCanvas *c0 = new TCanvas("c0","Canvas 0",400,10,600,700);
456 fPtpi->GetHistogram()->SetXTitle("p_{T} (GeV)");
459 fPtka->GetHistogram()->SetXTitle("p_{T} (GeV)");