]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/totEt/AliAnalysisEtCommon.cxx
Operating mode changes: Regular/Light/SuperLight Cascade output storage for Pb-Pb
[u/mrichter/AliRoot.git] / PWGLF / totEt / AliAnalysisEtCommon.cxx
1 //Created by Christine Nattrass, Rebecca Scott, Irakli Martashvili
2 //University of Tennessee at Knoxville
3 //
4 // This class is designed for the analysis of the hadronic component of 
5 // transverse energy.  It is used by AliAnalysisTaskHadEt.
6 // This gets information about the hadronic component of the transverse energy 
7 // from tracks reconstructed in an event
8 // it has daughters, AliAnalysisEtCommonMonteCarlo and 
9 // AliAnalysisEtCommonReconstructed which loop over either Monte Carlo data or 
10 // real data to get Et
11
12 #include "AliAnalysisEtCommon.h"
13 #include "TMath.h"
14 #include "TList.h"
15 #include "TH1F.h"
16 #include "TH2F.h"
17 #include "TF1.h"
18 #include <iostream>
19 #include "AliAnalysisEtCuts.h"
20 #include "AliMCEvent.h"
21 #include "AliVEvent.h"
22 #include "AliStack.h"
23 #include "AliESDtrackCuts.h"
24 #include "TDatabasePDG.h"
25 #include "TParticle.h"
26 #include "Rtypes.h"
27 #include "AliPDG.h"
28
29 using namespace std;
30
31 ClassImp(AliAnalysisEtCommon);
32 //These are from the PDG database but by making them static the code is a bit more efficient and has no problems running with the plugin
33 //Cuts are defined in $ROOTSYS/etc/pdg_table.txt
34
35 Float_t AliAnalysisEtCommon::fgPionMass = 0.13957;
36 Float_t AliAnalysisEtCommon::fgKaonMass = 0.493677;
37 Float_t AliAnalysisEtCommon::fgProtonMass = 0.938272;
38 Float_t AliAnalysisEtCommon::fgElectronMass = 0.000510999;
39 Int_t AliAnalysisEtCommon::fgPiPlusCode = 211;
40 Int_t AliAnalysisEtCommon::fgPiMinusCode = -211;
41 Int_t AliAnalysisEtCommon::fgKPlusCode = 321;
42 Int_t AliAnalysisEtCommon::fgKMinusCode = -321;
43 Int_t AliAnalysisEtCommon::fgProtonCode = 2212;
44 Int_t AliAnalysisEtCommon::fgAntiProtonCode = -2212;
45 Int_t AliAnalysisEtCommon::fgLambdaCode = 3122;
46 Int_t AliAnalysisEtCommon::fgAntiLambdaCode = -3122;
47 Int_t AliAnalysisEtCommon::fgK0Code = 311;
48 Int_t AliAnalysisEtCommon::fgK0SCode = 310;
49 Int_t AliAnalysisEtCommon::fgOmegaCode = 3334;
50 Int_t AliAnalysisEtCommon::fgAntiOmegaCode = -3334;
51 Int_t AliAnalysisEtCommon::fgXi0Code = 3322;
52 Int_t AliAnalysisEtCommon::fgAntiXi0Code = -3322;
53 Int_t AliAnalysisEtCommon::fgXiCode = 3312;
54 Int_t AliAnalysisEtCommon::fgAntiXiCode = -3312;
55 Int_t AliAnalysisEtCommon::fgSigmaCode = 3112;
56 Int_t AliAnalysisEtCommon::fgAntiSigmaCode = -3112;
57 Int_t AliAnalysisEtCommon::fgK0LCode = 130;
58 Int_t AliAnalysisEtCommon::fgNeutronCode = 2112;
59 Int_t AliAnalysisEtCommon::fgAntiNeutronCode = -2112;
60 Int_t AliAnalysisEtCommon::fgEPlusCode = -11;
61 Int_t AliAnalysisEtCommon::fgEMinusCode = 11;
62 Int_t AliAnalysisEtCommon::fgMuPlusCode = -13;
63 Int_t AliAnalysisEtCommon::fgMuMinusCode = 13;
64 Int_t AliAnalysisEtCommon::fgGammaCode = 22;
65 Int_t AliAnalysisEtCommon::fgPi0Code = 111;
66 Int_t AliAnalysisEtCommon::fgEtaCode = 221;
67 Int_t AliAnalysisEtCommon::fgOmega0Code = 223;
68
69
70
71 Float_t AliAnalysisEtCommon::fgPtTPCCutOff = 0.15;
72 Float_t AliAnalysisEtCommon::fgPtITSCutOff = 0.10;
73
74 AliAnalysisEtCommon::AliAnalysisEtCommon() : TObject()
75                                            ,fHistogramNameSuffix("")
76                                            ,fCuts(0)
77                                            ,fDataSet(2010)
78                                            ,fV0ScaleDataSet(2010)
79                                            ,fEsdtrackCutsITSTPC(0)
80                                            ,fEsdtrackCutsTPC(0)
81                                            ,fEsdtrackCutsITS(0)
82                                            ,fK0PythiaD6T(0)
83                                            ,fLambdaPythiaD6T(0)
84                                            ,fAntiLambdaPythiaD6T(0)
85                                            ,fK0Data(0)
86                                            ,fLambdaData(0)
87                                            ,fAntiLambdaData(0)
88                                            ,fLambdaEnhancement(0)
89                                            ,fProtonEnhancement(0)
90                                            ,fCentralityMethod("V0M")
91                                            ,fNCentBins(21)
92                                            ,fCentBin(-1)
93 {//default constructor
94
95 }
96
97 AliAnalysisEtCommon::~AliAnalysisEtCommon()
98 {//destructor
99   delete fCuts;
100   delete fEsdtrackCutsITSTPC;
101   delete fEsdtrackCutsITS;
102   delete fEsdtrackCutsTPC;
103   delete fK0PythiaD6T;
104   delete fLambdaPythiaD6T;
105   delete fAntiLambdaPythiaD6T;
106   delete fK0Data;
107   delete fLambdaData;
108   delete fAntiLambdaData;
109   delete fLambdaEnhancement;
110   delete fProtonEnhancement;
111 }
112
113 Int_t AliAnalysisEtCommon::AnalyseEvent(AliVEvent */*event*/)
114 { //this line is basically here to eliminate a compiler warning that event is not used.  Making it a virtual function did not work with the plugin.
115 //  cout << "This event has " << event->GetNumberOfTracks() << " tracks" << endl;
116   ResetEventValues();
117   return 0;
118 }
119
120
121 void AliAnalysisEtCommon::Init()
122 {// clear variables, set up cuts and PDG info
123   // LevyPt function described in LevyFitEvaluate below
124   //parameter 0 = dNdy
125   //parameter 1 = temp
126   //parameter 2 = power
127   if(fK0PythiaD6T) delete fK0PythiaD6T;
128   if(fLambdaPythiaD6T) delete fLambdaPythiaD6T;
129   if(fAntiLambdaPythiaD6T) delete fAntiLambdaPythiaD6T;
130   if(fK0Data) delete fK0Data;
131   if(fLambdaData) delete fLambdaData;
132   if(fAntiLambdaData) delete fAntiLambdaData;
133
134   fK0PythiaD6T = new TF1("K0PythiaD6T", LevyPtEvaluate, 0,50,4);
135   fLambdaPythiaD6T = new TF1("LambdaPythiaD6T", LevyPtEvaluate,0,50,4);
136   fAntiLambdaPythiaD6T = new TF1("AntiLambdaPythiaD6T", LevyPtEvaluate,0,50,4);
137   fK0Data = new TF1("K0Data", LevyPtEvaluate,0,50,4);
138   fLambdaData = new TF1("LambdaData", LevyPtEvaluate,0,50,4);
139   fAntiLambdaData = new TF1("AntiLambdaData", LevyPtEvaluate,0,50,4);
140
141   fK0PythiaD6T->FixParameter(3,0.493677);
142   fK0Data->FixParameter(3,0.493677);
143   fLambdaPythiaD6T->FixParameter(3,1.115683);
144   fAntiLambdaPythiaD6T->FixParameter(3,1.115683);
145   fLambdaData->FixParameter(3,1.115683);
146   fAntiLambdaData->FixParameter(3,1.115683);
147   if(fV0ScaleDataSet==2009){
148     //These data are from the ALICE 900 GeV p+p paper
149     //from email from Boris 10/31/10
150 //                              dN/dy                   Tlevy                   n                       <pt>     extra(%)       chi2/ndf
151 // K0S
152 // D6T-109      0.1437±0.0002    0.1497±0.0004   6.94±0.04       0.5765  14.24   162.5/24=6.8
153 // CSC-306      0.1938±0.0002    0.1363±0.0003   8.84±0.06       0.4950  16.80   1075.1/24=44.8
154 // P0-320               0.1395±0.0002    0.1440±0.0004   6.43±0.04       0.5778  14.56   387.0/24=16.1
155 // PHOJET               0.1967±0.0002    0.1292±0.0003   7.01±0.04       0.5115  16.90   1735.4/24=72.3
156
157 // Lambda
158 // D6T-109      0.0213±0.0001    0.1315±0.0011   4.60±0.05       0.7889  48.25   200.1/26=7.7
159 // CSC-306      0.0279±0.0001    0.1074±0.0007   5.71±0.07       0.6131  59.25   71.8/26=2.8
160 // P0-320               0.0161±0.0001    0.1125±0.0010   4.53±0.05       0.7126  53.61   171.1/26=6.6
161 // PHOJET               0.0190±0.0001    0.1013±0.0008   4.70±0.05       0.6469  58.14   162.2/26=6.2
162     //paper available at http://arxiv.org/pdf/1012.3257v2
163     //so that if I'm working at home and don't have library access I can read the paper...
164     //While all of these functions are labeled D6T, they are actually Perugia-0.  I did this to avoid the complication of renaming all of the variables.
165     fK0PythiaD6T->SetParameter(0,0.1395);
166     fK0PythiaD6T->SetParameter(1,0.1440);
167     fK0PythiaD6T->SetParameter(2,6.43);
168     fLambdaPythiaD6T->SetParameter(0,0.0161);
169     fLambdaPythiaD6T->SetParameter(1,0.1125);
170     fLambdaPythiaD6T->SetParameter(2,4.53);
171     fAntiLambdaPythiaD6T->SetParameter(0,0.0161);
172     fAntiLambdaPythiaD6T->SetParameter(1,0.1125);
173     fAntiLambdaPythiaD6T->SetParameter(2,4.53);
174     fK0Data->SetParameter(0,0.184);
175     fK0Data->SetParameter(1,0.168);
176     fK0Data->SetParameter(2,6.6);
177     fLambdaData->SetParameter(0,0.048);
178     fLambdaData->SetParameter(1,0.229);
179     fLambdaData->SetParameter(2,10.8);
180     fAntiLambdaData->SetParameter(0,0.047);
181     fAntiLambdaData->SetParameter(1,0.210);
182     fAntiLambdaData->SetParameter(2,9.2);
183   }
184   else{
185   //if(fV0ScaleDataSet==2010 ||fV0ScaleDataSet==20100 ||fV0ScaleDataSet==20111||fV0ScaleDataSet==2013){
186     //These data are from the CMS analysis note on 7 TeV spectra
187     //http://cdsweb.cern.ch/record/1279344/files/QCD-10-007-pas.pdf
188     //Note the CMS parameterization of the Levy function differs from the ALICE parameterization by a constant.
189     //CMS does not list the overall constant in their fit, the ratios of the dN/dy(y=0) is used.
190     //While all of these functions are labeled D6T, they are actually Perugia-0.  I did this to avoid the complication of renaming all of the variables.
191     fK0PythiaD6T->SetParameter(0,0.72);//K0S in CMS analysis note Fig. 10 - Perugia 0 and D6T level comparable
192     fK0PythiaD6T->SetParameter(1,0.168);
193     fK0PythiaD6T->SetParameter(2,5.39);
194     fLambdaPythiaD6T->SetParameter(0,0.54/1.12);//from looking at the difference in lambdas in Perugia 0 and D6T in fig. 10
195     fLambdaPythiaD6T->SetParameter(1,0.163);
196     fLambdaPythiaD6T->SetParameter(2,4.64);
197     fAntiLambdaPythiaD6T->SetParameter(0,0.54/1.12);//
198     fAntiLambdaPythiaD6T->SetParameter(1,0.163);
199     fAntiLambdaPythiaD6T->SetParameter(2,4.64);
200 //     fK0PythiaD6T->SetParameter(0,0.72);
201 //     fK0PythiaD6T->SetParameter(1,0.183);
202 //     fK0PythiaD6T->SetParameter(2,5.71);
203 //     fLambdaPythiaD6T->SetParameter(0,0.54);
204 //     fLambdaPythiaD6T->SetParameter(1,0.216);
205 //     fLambdaPythiaD6T->SetParameter(2,5.11);
206 //     fAntiLambdaPythiaD6T->SetParameter(0,0.54);
207 //     fAntiLambdaPythiaD6T->SetParameter(1,0.216);
208 //     fAntiLambdaPythiaD6T->SetParameter(2,5.11);
209     fK0Data->SetParameter(0,1.0);
210     fK0Data->SetParameter(1,0.215);
211     fK0Data->SetParameter(2,6.79);
212     fLambdaData->SetParameter(0,1.0);
213     fLambdaData->SetParameter(1,0.290);
214     fLambdaData->SetParameter(2,9.28);
215     fAntiLambdaData->SetParameter(0,1.0);
216     fAntiLambdaData->SetParameter(1,0.290);
217     fAntiLambdaData->SetParameter(2,9.28);
218   }
219   if(fLambdaEnhancement) delete fLambdaEnhancement;
220   fLambdaEnhancement = new TF1("fLambdaEnhancement","([0]*pow(x,[1])*exp(-pow(x/[2],[3])))/([4]*exp(-pow([5]/x,[6]))+[7]*x)",0,50);
221    fLambdaEnhancement->SetParameter(0,0.5630487);
222    fLambdaEnhancement->SetParameter(1,1.388818);
223    fLambdaEnhancement->SetParameter(2,3.954147);
224    fLambdaEnhancement->SetParameter(3,3.443772);
225    fLambdaEnhancement->SetParameter(4,2.844288);
226    fLambdaEnhancement->SetParameter(5,2);
227    fLambdaEnhancement->SetParameter(6,0.4747893);
228    fLambdaEnhancement->SetParameter(7,-0.2250856);
229    if(fProtonEnhancement) delete fProtonEnhancement;
230    fProtonEnhancement = new TF1("fProtonEnhancement","[0]*pow(x,[1])*exp(-pow(x/[2],[3]))/([4]+[5]*x)",0,50);
231    fProtonEnhancement->SetParameter(0,0.5630487*1.6);
232    fProtonEnhancement->SetParameter(1,1.388818);
233    fProtonEnhancement->SetParameter(2,3.954147/1.5);
234    fProtonEnhancement->SetParameter(3,3.443772/2.5);
235    fProtonEnhancement->SetParameter(4,0.5);
236    fProtonEnhancement->SetParameter(5,-.03);
237 }
238
239 void AliAnalysisEtCommon::ResetEventValues()
240 {//Resets event values of et to zero
241   
242   if (!fCuts) { // some Init's needed
243     //cout << __FILE__ << ":" << __LINE__ << " : Init " << endl;
244     if (!fCuts) {
245       cout << " setting up Cuts " << endl;
246       fCuts = new AliAnalysisEtCuts();
247     }
248   }
249 }
250
251 Double_t AliAnalysisEtCommon::LevyPtEvaluate(const Double_t *pt, 
252                                              const Double_t *par) 
253 {//LevyPt function for TF1's
254   
255   Double_t lMass  = par[3];
256   Double_t ldNdy  = par[0];
257   Double_t l2pi   = 2*TMath::Pi();
258   Double_t lTemp = par[1];
259   Double_t lPower = par[2];
260   
261   Double_t lBigCoef = ((lPower-1)*(lPower-2)) / (l2pi*lPower*lTemp*(lPower*lTemp+lMass*(lPower-2)));
262   Double_t lInPower = 1 + (TMath::Sqrt(pt[0]*pt[0]+lMass*lMass)-lMass) / (lPower*lTemp);
263   
264   return ldNdy * pt[0] * lBigCoef * TMath::Power(lInPower,(-1)*lPower);
265 }
266
267
268 Float_t AliAnalysisEtCommon::Et(TParticle *part, float mass){//function to calculate et in the same way as it would be calculated in a calorimeter
269   if(mass+1000<0.01){//if no mass given return default.  The default argument is -1000
270     if(TMath::Abs(part->GetPDG(0)->PdgCode())==2212 || TMath::Abs(part->GetPDG(0)->PdgCode())==2112){
271       if(part->GetPDG(0)->PdgCode()==-2212 || part->GetPDG(0)->PdgCode()==-2112){//antiproton or antineutron
272         //for antinucleons we specifically want to return the kinetic energy plus twice the rest mass
273         return (part->Energy()+part->GetMass())*TMath::Sin(part->Theta());
274       }
275       if(part->GetPDG(0)->PdgCode()==2212 || part->GetPDG(0)->PdgCode()==2112){//proton or neutron
276         //for nucleons we specifically want to return the kinetic energy only
277         return (part->Energy()-part->GetMass())*TMath::Sin(part->Theta());
278       }
279     }
280     else{//otherwise go to the default
281       return part->Energy()*TMath::Sin(part->Theta());
282     }
283   }
284   else{//otherwise use the mass that was given
285     return (TMath::Sqrt(TMath::Power(part->P(),2.0)+TMath::Power(mass,2.0)))*TMath::Sin(part->Theta());
286   }
287   return 0.0;
288 }
289 Float_t AliAnalysisEtCommon::Et(Float_t p, Float_t theta, Int_t pid, Short_t charge) const {//function to calculate et in the same way as it would be calculated in a calorimeter
290   if(pid==fgPiPlusCode || pid==fgPiMinusCode){//Nothing special for pions
291     return TMath::Sqrt(p*p + fgPionMass*fgPionMass) * TMath::Sin(theta);
292   }
293   if(pid==fgKPlusCode || pid==fgKMinusCode){//Nothing special for kaons
294     return TMath::Sqrt(p*p + fgKaonMass*fgKaonMass) * TMath::Sin(theta);
295   }
296   if(pid==fgEPlusCode || pid==fgEMinusCode){//Nothing special for electrons
297     return TMath::Sqrt(p*p + fgElectronMass*fgElectronMass) * TMath::Sin(theta);
298   }
299   if(pid==fgProtonCode || pid==fgAntiProtonCode){//But for protons we must be careful...
300     if(charge<0.0){//antiprotns: kinetic energy plus twice the rest mass
301       return (TMath::Sqrt(p*p + fgProtonMass*fgProtonMass) + fgProtonMass) * TMath::Sin(theta);
302     }
303     if(charge>0.0){//antiprotns: kinetic energy only
304       return (TMath::Sqrt(p*p + fgProtonMass*fgProtonMass) - fgProtonMass) * TMath::Sin(theta);
305     }
306   }
307   cerr<<"Uh-oh!  Et not set properly!"<<endl;
308   return 0.0;
309 }
310
311 Float_t AliAnalysisEtCommon::K0Weight(Float_t pt){
312   Float_t data = fK0Data->Eval(pt);
313   Float_t mc = fK0PythiaD6T->Eval(pt);
314   return data/mc;
315 }
316 Float_t AliAnalysisEtCommon::LambdaWeight(Float_t pt){
317   Float_t data = fLambdaData->Eval(pt);
318   Float_t mc = fLambdaPythiaD6T->Eval(pt);
319   return data/mc;
320 }
321 Float_t AliAnalysisEtCommon::AntiLambdaWeight(Float_t pt){
322   Float_t data = fAntiLambdaData->Eval(pt);
323   Float_t mc = fAntiLambdaPythiaD6T->Eval(pt);
324   return data/mc;
325 }
326
327 Float_t AliAnalysisEtCommon::LambdaBaryonEnhancement(Float_t pt){
328   if(pt<0.8) return 1.0;
329   return fLambdaEnhancement->Eval(pt);
330 };//Function which gives the factor to reweigh a lambda or antilambda so it roughly matches baryon enhancement seen at RHIC
331 Float_t AliAnalysisEtCommon::ProtonBaryonEnhancement(Float_t pt){
332   if(pt<0.8) return 1.0;
333   return fProtonEnhancement->Eval(pt);
334 }//Function which gives the factor to reweigh a lambda or antilambda so it roughly matches baryon enhancement seen at RHIC