]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EVGEN/AliGenHIJINGpara.cxx
Removing leftover return
[u/mrichter/AliRoot.git] / EVGEN / AliGenHIJINGpara.cxx
CommitLineData
790bbabf 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
88cb7938 16/* $Id$ */
675e9664 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
790bbabf 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
971816d4 40#include <TArrayF.h>
ac3faee4 41#include <TCanvas.h>
2067f62c 42#include <TClonesArray.h>
116cbefd 43#include <TDatabasePDG.h>
44#include <TF1.h>
cc692c80 45#include <TH1.h>
116cbefd 46#include <TPDGCode.h>
ac3faee4 47#include <TParticle.h>
828fff0d 48#include <TROOT.h>
ac3faee4 49#include <TVirtualMC.h>
116cbefd 50
51#include "AliConst.h"
52#include "AliDecayer.h"
116cbefd 53#include "AliGenEventHeader.h"
54#include "AliGenHIJINGpara.h"
e087ab36 55#include "AliLog.h"
116cbefd 56#include "AliRun.h"
971816d4 57
790bbabf 58ClassImp(AliGenHIJINGpara)
59
198bb1c7 60
790bbabf 61
62//_____________________________________________________________________________
4a33c50d 63static Double_t ptpi(const Double_t *px, const Double_t *)
790bbabf 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 //
cc692c80 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.;
790bbabf 76 Double_t y, y1, xmpi2, ynorm, a;
cc692c80 77 Double_t x = *px;
790bbabf 78 //
cc692c80 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;
790bbabf 83 if (x > kxlim)
cc692c80 84 y = a * TMath::Power(kp0 / (kp0 + x), kxn);
790bbabf 85 else
cc692c80 86 y = kb* TMath::Exp(-sqrt(x * x + xmpi2) / kt);
87
790bbabf 88 return y*x;
89}
90
91//_____________________________________________________________________________
92static 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//_____________________________________________________________________________
112static 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//_____________________________________________________________________________
122static 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//_____________________________________________________________________________
141static 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//_____________________________________________________________________________
160AliGenHIJINGpara::AliGenHIJINGpara()
1c56e311 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)
790bbabf 172{
173 //
174 // Default constructor
175 //
971816d4 176 SetCutVertexZ();
6b236b6f 177 SetPtRange();
790bbabf 178}
179
180//_____________________________________________________________________________
181AliGenHIJINGpara::AliGenHIJINGpara(Int_t npart)
1c56e311 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)
790bbabf 193{
194 //
195 // Standard constructor
196 //
edbba5db 197 fName="HIJINGpara";
790bbabf 198 fTitle="HIJING Parametrisation Particle Generator";
971816d4 199 SetCutVertexZ();
6b236b6f 200 SetPtRange();
1c56e311 201}
202
790bbabf 203//_____________________________________________________________________________
204AliGenHIJINGpara::~AliGenHIJINGpara()
205{
206 //
207 // Standard destructor
208 //
209 delete fPtpi;
210 delete fPtka;
211 delete fETApic;
212 delete fETAkac;
213}
214
215//_____________________________________________________________________________
216void 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)));
f25521b2 225 fPtpi = new TF1("ptpi",&ptpi,0,20,0);
828fff0d 226 gROOT->GetListOfFunctions()->Remove(fPtpi);
f25521b2 227 fPtka = new TF1("ptka",&ptka,0,20,0);
828fff0d 228 gROOT->GetListOfFunctions()->Remove(fPtka);
cc692c80 229 fPtpi->SetNpx(1000);
230 fPtka->SetNpx(1000);
790bbabf 231 fETApic = new TF1("etapic",&etapic,etaMin,etaMax,0);
828fff0d 232 gROOT->GetListOfFunctions()->Remove(fETApic);
790bbabf 233 fETAkac = new TF1("etakac",&etakac,etaMin,etaMax,0);
828fff0d 234 gROOT->GetListOfFunctions()->Remove(fETAkac);
f25521b2 235
bbbc190a 236 TF1 etaPic0("etaPic0",&etapic,-7,7,0);
237 TF1 etaKac0("etaKac0",&etakac,-7,7,0);
f25521b2 238
bbbc190a 239 TF1 ptPic0("ptPic0",&ptpi,0.,15.,0);
240 TF1 ptKac0("ptKac0",&ptka,0.,15.,0);
f25521b2 241
aecf26a3 242 Float_t intETApi = etaPic0.Integral(-0.5, 0.5);
243 Float_t intETAka = etaKac0.Integral(-0.5, 0.5);
f25521b2 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
aecf26a3 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));
f25521b2 252 Float_t ptFrac = intPtSel/intPt;
253
254// Fraction of events corresponding to the selected eta-range
aecf26a3 255 Float_t intETASel = (scalePi*etaPic0.Integral(etaMin, etaMax)+
256 scaleKa*etaKac0.Integral(etaMin, etaMax));
f25521b2 257// Fraction of events corresponding to the selected phi-range
258 Float_t phiFrac = (fPhiMax-fPhiMin)/2/TMath::Pi();
259
0af12c00 260
647f0914 261 fParentWeight = (intETASel*ptFrac*phiFrac) / Float_t(fNpart);
790bbabf 262
0af12c00 263 if (fAnalog != 0) {
264 fPtWgtPi = (fPtMax - fPtMin) / fPtpi->Integral(0., 20.);
265 fPtWgtKa = (fPtMax - fPtMin) / fPtka->Integral(0., 20.);
647f0914 266 fParentWeight = (intETASel*phiFrac) / Float_t(fNpart);
0af12c00 267 }
268
269
e087ab36 270 AliInfo(Form("The number of particles in the selected kinematic region corresponds to %f percent of a full event",
647f0914 271 100./ fParentWeight));
5b0e4c0f 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) {
e087ab36 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));
5b0e4c0f 279 }
2067f62c 280//
281//
2942f542 282 if (fPi0Decays && TVirtualMC::GetMC())
283 fDecayer = TVirtualMC::GetMC()->GetDecayer();
ebf41d7c 284
285 if (fPi0Decays)
286 {
287 fDecayer->SetForceDecay(kNeutralPion);
288 fDecayer->Init();
289 }
290
790bbabf 291}
292
2067f62c 293
790bbabf 294//_____________________________________________________________________________
295void AliGenHIJINGpara::Generate()
296{
297 //
298 // Generate one trigger
299 //
300
301
302 const Float_t kRaKpic=0.14;
303 const Float_t kBorne=1/(1+kRaKpic);
304 Float_t polar[3]= {0,0,0};
305 //
306 const Int_t kPions[3] = {kPi0, kPiPlus, kPiMinus};
307 const Int_t kKaons[4] = {kK0Long, kK0Short, kKPlus, kKMinus};
308 //
309 Float_t origin[3];
21391258 310 Float_t time;
0af12c00 311 Float_t pt, pl, ptot, wgt;
790bbabf 312 Float_t phi, theta;
313 Float_t p[3];
1e66857c 314 Int_t i, part, j;
790bbabf 315 //
316 TF1 *ptf;
317 TF1 *etaf;
318 //
319 Float_t random[6];
320 //
321 for (j=0;j<3;j++) origin[j]=fOrigin[j];
21391258 322 time = fTimeOrigin;
971816d4 323
324 if(fVertexSmear == kPerEvent) {
c8f7f6f9 325 Vertex();
326 for (j=0; j < 3; j++) origin[j] = fVertex[j];
21391258 327 time = fTime;
971816d4 328 } // if kPerEvent
329 TArrayF eventVertex;
330 eventVertex.Set(3);
331 eventVertex[0] = origin[0];
332 eventVertex[1] = origin[1];
333 eventVertex[2] = origin[2];
21391258 334 Float_t eventTime = time;
0af12c00 335
790bbabf 336 for(i=0;i<fNpart;i++) {
337 while(1) {
0af12c00 338 Rndm(random,4);
790bbabf 339 if(random[0]<kBorne) {
340 part=kPions[Int_t (random[1]*3)];
341 ptf=fPtpi;
2067f62c 342 etaf=fETApic;
0af12c00 343 wgt = fPtWgtPi;
790bbabf 344 } else {
345 part=kKaons[Int_t (random[1]*4)];
346 ptf=fPtka;
347 etaf=fETAkac;
0af12c00 348 wgt = fPtWgtKa;
790bbabf 349 }
350 phi=fPhiMin+random[2]*(fPhiMax-fPhiMin);
351 theta=2*TMath::ATan(TMath::Exp(-etaf->GetRandom()));
352 if(theta<fThetaMin || theta>fThetaMax) continue;
0af12c00 353
354 if (fAnalog == 0) {
355 pt = ptf->GetRandom();
356 } else {
357 pt = fPtMin + random[3] * (fPtMax - fPtMin);
358 }
359
360
790bbabf 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);
366 p[2]=pl;
aee8290b 367 if(fVertexSmear==kPerTrack) {
65fb704d 368 Rndm(random,6);
790bbabf 369 for (j=0;j<3;j++) {
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]));
372 }
21391258 373 Rndm(random,2);
374 time = fTimeOrigin + fOsigma[2]/TMath::Ccgs()*
375 TMath::Cos(2*random[0]*TMath::Pi())*
376 TMath::Sqrt(-2*TMath::Log(random[1]));
790bbabf 377 }
a142c195 378
0af12c00 379 if (fAnalog == 0) {
380 wgt = fParentWeight;
381 } else {
382 wgt *= (fParentWeight * ptf->Eval(pt));
383 }
384
385
1ef5f15d 386 if (part == kPi0 && fPi0Decays){
2067f62c 387//
388// Decay pi0 if requested
21391258 389 PushTrack(0,-1,part,p,origin,polar,time,kPPrimary,fNt,wgt);
1e66857c 390 KeepTrack(fNt);
21391258 391 DecayPi0(origin, p, time);
1ef5f15d 392 } else {
a142c195 393 // printf("fNt %d", fNt);
21391258 394 PushTrack(fTrackIt,-1,part,p,origin,polar,time,kPPrimary,fNt,wgt);
a142c195 395
1e66857c 396 KeepTrack(fNt);
1ef5f15d 397 }
398
790bbabf 399 break;
400 }
1e66857c 401 SetHighWaterMark(fNt);
790bbabf 402 }
1e66857c 403//
404
971816d4 405// Header
406 AliGenEventHeader* header = new AliGenEventHeader("HIJINGparam");
407// Event Vertex
408 header->SetPrimaryVertex(eventVertex);
21391258 409 header->SetInteractionTime(eventTime);
a142c195 410 header->SetNProduced(fNpartProd);
45c05df7 411 if (fContainer) {
412 header->SetName(fName);
413 fContainer->AddHeader(header);
414 } else {
415 gAlice->SetGenEventHeader(header);
416 }
790bbabf 417}
418
6b236b6f 419void AliGenHIJINGpara::SetPtRange(Float_t ptmin, Float_t ptmax) {
2067f62c 420 AliGenerator::SetPtRange(ptmin, ptmax);
6b236b6f 421}
790bbabf 422
21391258 423void AliGenHIJINGpara::DecayPi0(Float_t* orig, Float_t * p, Float_t time)
2067f62c 424{
425//
426// Decay the pi0
427// and put decay products on the stack
428//
429 static TClonesArray *particles;
430 if(!particles) particles = new TClonesArray("TParticle",1000);
431//
432 const Float_t kMass = TDatabasePDG::Instance()->GetParticle(kPi0)->Mass();
56fcdfe6 433 Float_t e = TMath::Sqrt(p[0] * p[0] + p[1] * p[1] + p[2] * p[2]+ kMass * kMass);
2067f62c 434//
435// Decay the pi0
436 TLorentzVector pmom(p[0], p[1], p[2], e);
437 fDecayer->Decay(kPi0, &pmom);
56fcdfe6 438
2067f62c 439//
440// Put decay particles on the stack
441//
2067f62c 442 Float_t polar[3] = {0., 0., 0.};
443 Int_t np = fDecayer->ImportParticles(particles);
a142c195 444 fNpartProd += (np-1);
8f64c44f 445 Int_t nt = 0;
2067f62c 446 for (Int_t i = 1; i < np; i++)
447 {
7b33c00a 448 TParticle* iParticle = (TParticle *) particles->At(i);
449 p[0] = iParticle->Px();
450 p[1] = iParticle->Py();
451 p[2] = iParticle->Pz();
2067f62c 452 Int_t part = iParticle->GetPdgCode();
1e66857c 453
21391258 454 PushTrack(fTrackIt, fNt, part, p, orig, polar, time, kPDecay, nt, fParentWeight);
1e66857c 455 KeepTrack(nt);
2067f62c 456 }
1e66857c 457 fNt = nt;
2067f62c 458}
198bb1c7 459
cc692c80 460void AliGenHIJINGpara::Draw( const char * /*opt*/)
461{
462 //
463 // Draw the pT and y Distributions
464 //
465 TCanvas *c0 = new TCanvas("c0","Canvas 0",400,10,600,700);
466 c0->Divide(2,1);
467 c0->cd(1);
468 fPtpi->Draw();
469 fPtpi->GetHistogram()->SetXTitle("p_{T} (GeV)");
470 c0->cd(2);
471 fPtka->Draw();
472 fPtka->GetHistogram()->SetXTitle("p_{T} (GeV)");
473
474}