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