]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EVGEN/AliGenTunedOnPbPb.cxx
Correction
[u/mrichter/AliRoot.git] / EVGEN / AliGenTunedOnPbPb.cxx
CommitLineData
995e33f0 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 *
995e33f0 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: AliGenTunedOnPbPb.cxx 51126 2013-08-19 13:37:49Z fnoferin $ */
17
18// Parameterisation based on 5.5 ATeV PbPb data
19// pi,K, p , K0, lambda, phi, Xi, Omega spectra, v2, v3 (no jets!)
20// Author: fnoferin@cern.ch
21
22#include <TArrayF.h>
23#include <TCanvas.h>
24#include <TClonesArray.h>
25#include <TDatabasePDG.h>
26#include <TF1.h>
27#include <TH1.h>
28#include <TPDGCode.h>
29#include <TParticle.h>
30#include <TROOT.h>
31#include <TVirtualMC.h>
32
33#include "AliConst.h"
34#include "AliDecayer.h"
35#include "AliGenEventHeaderTunedPbPb.h"
36#include "AliLog.h"
37#include "AliRun.h"
38#include "AliGenTunedOnPbPb.h"
39
40ClassImp(AliGenTunedOnPbPb)
41
42// set default parameters for 10-20% centrality
43Int_t AliGenTunedOnPbPb::fgPdgInput[fgNspecies] = {211,-211,111,321,-321,2212,-2212,310,3122,-3122,333,3312,-3312,3334,-3334};
44Float_t AliGenTunedOnPbPb::fgMult[fgNspecies] = {450,450,450,70,70,21,21,70,20,20,8,2.4,2.4,0.4,0.4};
45
46Float_t AliGenTunedOnPbPb::fgV3Overv2 = 6.25000000000000000e-01;
47Float_t AliGenTunedOnPbPb::fgEventplane=0;
48
49TF1 *AliGenTunedOnPbPb::fgV2 = NULL;
50
51//_____________________________________________________________________________
52AliGenTunedOnPbPb::AliGenTunedOnPbPb()
53 :AliGenerator(),
54 fCmin(0.),
55 fCmax(100.),
56 fChangeWithCentrality(kFALSE),
57 fYMaxValue(2.0)
58{
59 //
60 // Default constructor
61 //
62 SetCutVertexZ();
63 SetPtRange();
64
65 for(Int_t i=0;i < fgNspecies;i++){
66 fgHSpectrum[i] = NULL;
67 fgHv2[i] = NULL;
68 }
69}
70
71//_____________________________________________________________________________
72AliGenTunedOnPbPb::~AliGenTunedOnPbPb()
73{
74 //
75 // Standard destructor
76 //
77}
78
79//_____________________________________________________________________________
80void AliGenTunedOnPbPb::Init()
81{
82 //
83 // Initialise the generator
84 //
85
86 // define histos
87}
88
89
90//_____________________________________________________________________________
91void AliGenTunedOnPbPb::Generate()
92{
93 //
94 // Generate one trigger
95 //
96
97 Float_t avCentr = (fCmin+fCmax)*0.5;
98
99 Float_t centrality = avCentr;
100
101 if(fChangeWithCentrality) centrality = fCmin + gRandom->Rndm()*(fCmax-fCmin);
102
103 SetParameters(centrality);
104
105 TMCProcess statusPdg[fgNspecies] = {kPPrimary,kPPrimary,kPPrimary,kPPrimary,kPPrimary,kPPrimary,kPPrimary,kPPrimary,kPPrimary,kPPrimary,kPPrimary,kPPrimary,kPPrimary,kPPrimary,kPPrimary};
106
107 Float_t parV2scaling[3] = {1,0.202538,-0.00214468};
108
109 Float_t scaleV2 = 1.0;
110
111 TDatabasePDG *pdgD = TDatabasePDG::Instance();
112
113 if(fChangeWithCentrality){
114 parV2scaling[0] = 1. / (1 + parV2scaling[1]*avCentr + parV2scaling[2]*avCentr*avCentr); // normalize to average centrality
115
116 scaleV2 = parV2scaling[0]*(1 + parV2scaling[1]*centrality + parV2scaling[2]*centrality*centrality); // apply the trand of v2 w.r.t. centrality
117 }
118
119 fgV2->SetParameter(2,fgV3Overv2);
120
121 Float_t psi = gRandom->Rndm()*TMath::Pi();
122 fgEventplane = psi;
123 Float_t psi3 = gRandom->Rndm()*TMath::Pi()*2/3;
124 fgV2->SetParameter(1,psi);
125 fgV2->SetParameter(3,psi3);
126
127 Int_t npart = 0;
128
129 Float_t origin0[3];
130 for (Int_t j=0;j<3;j++) origin0[j]=fOrigin[j];
131 Float_t time;
132 time = fTimeOrigin;
133 if(fVertexSmear == kPerEvent) {
134 Vertex();
135 for (Int_t j=0; j < 3; j++) origin0[j] = fVertex[j];
136 time = fTime;
137 } // if kPerEvent
138
139 printf("Generate event with centrality = %3.1f%c, |y|<%4.1f\n",centrality,'%',fYMaxValue);
140
141 for(Int_t isp=0;isp < fgNspecies;isp++){
142 if(! fgHSpectrum[isp]) continue;
143
144 printf("Total number of %i = %f\n",fgPdgInput[isp],fgMult[isp]*2*fYMaxValue);
145
146 for(Int_t ipart =0; ipart < fgMult[isp]*2*fYMaxValue; ipart++){
147 Int_t pdg = fgPdgInput[isp];
148
149 Double_t y = gRandom->Rndm()*2*fYMaxValue - fYMaxValue;
150 Double_t pt = fgHSpectrum[isp]->GetRandom();
151 if(fgHv2[isp]) fgV2->SetParameter(0,fgHv2[isp]->Interpolate(pt) * scaleV2);
152 else fgV2->SetParameter(0,0.);
153 Double_t mass = pdgD->GetParticle(pdg)->Mass();
154 Double_t mt2 = pt*pt + mass*mass;
155 Double_t phi = fgV2->GetRandom(-TMath::Pi(),TMath::Pi());
156 Double_t px = pt*TMath::Cos(phi);
157 Double_t py = pt*TMath::Sin(phi);
158 y = TMath::TanH(y);
159 Double_t pz = y*TMath::Sqrt(mt2)/TMath::Sqrt(1-y*y);
160// Double_t etot = TMath::Sqrt(mt2 + pz*pz);
161 Float_t p[3] = {px,py,pz};
162 Float_t polar[3] = {0.,0.,0.};
163
164 PushTrack(1, -1, pdg, p, origin0, polar, time, statusPdg[isp], npart, 1., 1);
165 KeepTrack(npart);
166 npart++;
167 }
168 }
169
170 TArrayF eventVertex;
171 eventVertex.Set(3);
172 eventVertex[0] = origin0[0];
173 eventVertex[1] = origin0[1];
174 eventVertex[2] = origin0[2];
175
176// Header
177 AliGenEventHeaderTunedPbPb* header = new AliGenEventHeaderTunedPbPb("tunedOnPbPb");
178// Event Vertex
179 header->SetPrimaryVertex(eventVertex);
180 header->SetInteractionTime(time);
181 header->SetNProduced(npart);
182 header->SetCentrality(centrality);
183 header->SetPsi2(psi);
184 header->SetPsi3(psi3);
185 gAlice->SetGenEventHeader(header);
186}
187
188void AliGenTunedOnPbPb::SetPtRange(Float_t ptmin, Float_t ptmax) {
189 AliGenerator::SetPtRange(ptmin, ptmax);
190}
191//_____________________________________________________________________________
192TH1F *AliGenTunedOnPbPb::GetMultVsCentrality(Int_t species){
193 char title[100];
194 snprintf(title,100,"pdg = %i;centrality;dN/dy",fgPdgInput[species]);
195 TH1F *h = new TH1F("multVsCentr",title,100,0,100);
196
197 for(Int_t i=1;i<=100;i++){
198 Float_t x = i+0.5;
199 SetParameters(x);
200 h->SetBinContent(i,fgMult[species]);
201 }
202
203 return h;
204}
205//_____________________________________________________________________________
206void AliGenTunedOnPbPb::SetParameters(Float_t centrality){
207
208 if(!fgV2) fgV2 = new TF1("fv2Par","(1 + 2*[0]*cos(2*(x-[1])) + 2*[0]*[2]*cos(3*(x-[3])))",-TMath::Pi(),TMath::Pi());
209
210 Float_t fr[9] = {0.,0.,0.,0.,0.,0.,0.,0.,0.};
211
212 if(centrality < 7.5){
213 fr[0] = (7.5 - centrality)/5.;
214 fr[1] = (centrality-2.5)/5.;
215 }
216 else if(centrality < 15){
217 fr[1] = (15-centrality)/7.5;
218 fr[2] = (centrality-7.5)/7.5;
219 }
220 else if(centrality < 25){
221 fr[2] = (25-centrality)/10.;
222 fr[3] = (centrality-15)/10.;
223 }
224 else if(centrality < 35){
225 fr[3] = (35-centrality)/10.;
226 fr[4] = (centrality-25)/10.;
227 }
228 else if(centrality < 45){
229 fr[4] = (45-centrality)/10.;
230 fr[5] = (centrality-35)/10.;
231 }
232 else if(centrality < 55){
233 fr[5] = (55-centrality)/10.;
234 fr[6] = (centrality-45)/10.;
235 }
236 else if(centrality < 65){
237 fr[6] = (65-centrality)/10.;
238 fr[7] = (centrality-55)/10.;
239 }
240 else if(centrality < 75){
241 fr[7] = (75-centrality)/10.;
242 fr[8] = (centrality-65)/10.;
243 }
244 else{
245 fr[8] = 1.0;
246 }
247
248 // parameters as a function of centrality
249 Float_t multCent[9*fgNspecies] = {680.,680.,680.,110.,110.,34.,34.,110.,28.,28.,11.5,3.1,3.1,0.5,0.5,
250 560.,560.,560.,92.,92.,28.,28.,92.,24.,24.,9.,2.7,2.7,0.45,0.45,
251 450.,450.,450.,70.,70.,21.,21.,70.,20.,20.,8.,2.4,2.4,0.4,0.4,
252 302.,302.,302.,47.,47.,14.5,14.5,47.,14.,14.,5.5,1.5,1.5,0.2,0.2,
253 200.,200.,200.,30.,30.,9.6,9.6,30.,9.,9.,3.5,0.9,0.9,0.08,0.08,
254 120.,120.,120.,18.,18.,6.1,6.1,18.,5.1,5.1,2.2,0.6,0.6,0.055,0.055,
255 70.,70.,70.,10.4,10.4,3.6,3.6,10.4,2.6,2.6,1.4,0.36,0.36,0.035,0.035,
256 36.,36.,36.,5.25,5.25,2.,2.,5.25,1.5,1.5,0.5,0.02,0.02,0.015,0.015,
257 17.,17.,17.,2.3,2.3,0.9,0.9,2.3,0.6,0.6,0.16,0.006,0.006,0.005,0.005};
258
259 Float_t v3Overv2Cent[9] = {1.2,0.82,0.625,0.5,0.45,0.4,0.37,0.3,0.3};
260
261 // set parameters for current centrality
262 for(Int_t i=0;i < fgNspecies;i++){
263 fgMult[i] = 0;
264
265 for(Int_t j=0;j < 9;j++){
266 fgMult[i] += fr[j]*multCent[i+j*fgNspecies];
267 }
268 }
269
270 fgV3Overv2 = 0;
271 for(Int_t j=0;j < 9;j++)
272 fgV3Overv2 += fr[j]*v3Overv2Cent[j];
273
274 if(centrality > 80){
275 for(Int_t i=0;i < fgNspecies;i++)
276 fgMult[i] /= TMath::Log(centrality-77.);
277 }
278
279}