]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/totEt/AliAnalysisHadEt.cxx
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / PWGLF / totEt / AliAnalysisHadEt.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, AliAnalysisHadEtMonteCarlo and 
9 // AliAnalysisHadEtReconstructed which loop over either Monte Carlo data or 
10 // real data to get Et
11
12 #include "AliAnalysisHadEt.h"
13 #include "TMath.h"
14 #include "TList.h"
15 #include "TH1F.h"
16 #include "TH2F.h"
17 #include <iostream>
18 #include "AliAnalysisEtCuts.h"
19 #include "AliMCEvent.h"
20 #include "AliVEvent.h"
21 #include "AliStack.h"
22 #include "AliESDtrackCuts.h"
23 #include "TDatabasePDG.h"
24 #include "TParticle.h"
25 #include "Rtypes.h"
26 #include "AliPDG.h"
27 #include "AliPIDResponse.h"
28 #include "AliTPCPIDResponse.h" 
29 #include "AliInputEventHandler.h"
30 #include "AliAnalysisManager.h"
31 #include "AliLog.h"
32 #include "AliCentrality.h"
33
34 using namespace std;
35
36 ClassImp(AliAnalysisHadEt);
37
38
39 Int_t AliAnalysisHadEt::fgnumOfEtaBins = 16;
40 Float_t AliAnalysisHadEt::fgEtaAxis[17]={-0.78, -0.7, -0.58, -0.46, -0.34, -0.22, -0.12, -0.06, -0.0, 0.06, 0.12, 0.22, 0.34, 0.46, 0.58, 0.7, 0.78};
41 Int_t AliAnalysisHadEt::fgNumOfPtBins = 111;
42 Float_t AliAnalysisHadEt::fgPtAxis[117]=
43   {0.0,0.01,0.02,0.03,0.04, 0.05, 0.06,0.07,0.08,0.09, 0.10,0.11, .12,0.13, .14,0.15, .16,0.17, .18,0.19,
44    0.2, .22, .24, .26, .28, 0.30, 0.32, .34, .36, .38, 0.40, .42, .44, .46, .48,
45    0.5, .52, .54, .56, .58, 0.60, 0.62, .64, .66, .68, 0.70, .72, .74, .76, .78,
46    .80, .82, .84, .86, .88, 0.90, 0.92, .94, .96, .98, 1.00,1.05, 1.1,1.15, 1.2,
47   1.25, 1.3,1.35,1.40,1.45, 1.50, 1.55, 1.6,1.65, 1.7, 1.75, 1.8,1.85, 1.9,1.95,
48    2.0, 2.2, 2.4, 2.6, 2.8, 3.00, 3.20, 3.4, 3.6, 3.8, 4.00, 4.2, 4.4, 4.6, 4.8,
49    5.0, 5.5, 6.0, 6.5, 7.0, 7.50, 8.00, 8.5, 9.0, 9.5, 10.0,12.0,14.0,16.0,18.0,
50   20.0,25.0,30.0,35.0,40.0, 45.0, 50.0};
51 Float_t AliAnalysisHadEt::fgResAxis[81] = {-0.150,-0.140,-0.130,-0.120,-0.110,-0.100,-0.090,-0.080,-0.070,-0.060,
52                                            -0.050,-0.045,-0.040,-0.035,-0.030,-0.025,-0.024,-0.023,-0.022,-0.021,
53                                            -0.020,-0.019,-0.018,-0.017,-0.016,-0.015,-0.014,-0.013,-0.012,-0.011,
54                                            -0.010,-0.009,-0.008,-0.007,-0.006,-0.005,-0.004,-0.003,-0.002,-0.001,
55                                            -0.000,0.001,0.002,0.003,0.004,0.005,0.006,0.007,0.008,0.009,
56                                            0.010,0.011,0.012,0.013,0.014,0.015,0.016,0.017,0.018,0.019,
57                                            0.020,0.021,0.022,0.023,0.024,0.025,0.030,0.035,0.040,0.045,
58                                            0.050,0.060,0.070,0.080,0.090,0.100,0.110,0.120,0.130,0.140,
59                                            0.150,};
60 // Float_t AliAnalysisHadEt::fgResAxis[31] = {-0.15,-0.14,-0.13,-0.12,-0.11,-0.1,-0.09,-0.08,-0.07,-0.06,-0.05,-0.04,-0.03,-0.02,-0.01,0.0,0.01,0.02,0.03,0.04,0.05,0.06,0.07,0.08,0.09,0.1,0.11,0.12,0.13,0.14,0.15};
61 Int_t AliAnalysisHadEt::fgNumOfResBins = 80;
62
63
64 AliAnalysisHadEt::AliAnalysisHadEt() : AliAnalysisEtCommon()
65                                      ,fPIDResponse(0)
66                                      ,fSumEt(0)
67                                      ,fSumEtAcc(0)
68                                      ,fTotEt(0)
69                                      ,fTotEtAcc(0)
70                                      ,fTotNeutralEt(0)
71                                      ,fTotNeutralEtAcc(0)
72                                      ,fTotChargedEt(0)
73                                      ,fTotChargedEtAcc(0)
74                                      ,fMultiplicity(0)
75                                      ,fChargedMultiplicity(0)
76                                      ,fNeutralMultiplicity(0)
77                                      ,fhistoList(0)
78                                      ,fGoodEvent(0)
79 {//default constructor
80
81 }
82
83 AliAnalysisHadEt::~AliAnalysisHadEt()
84 {//destructor
85   delete fEsdtrackCutsITSTPC;
86   delete fEsdtrackCutsITS;
87   delete fEsdtrackCutsTPC;
88   delete fPIDResponse;
89 }
90
91 Int_t AliAnalysisHadEt::AnalyseEvent(AliVEvent *event)
92 { //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.
93   AliAnalysisEtCommon::AnalyseEvent(event);
94   ResetEventValues();
95   return 0;
96 }
97
98
99 void AliAnalysisHadEt::Init()
100 {// clear variables, set up cuts and PDG info
101   AliAnalysisEtCommon::Init();
102   
103 }
104
105 void AliAnalysisHadEt::ResetEventValues()
106 {//Resets event values of et to zero
107   AliAnalysisEtCommon::ResetEventValues();
108   fTotEt = 0;
109   fTotEtAcc = 0;
110   fTotNeutralEt = 0;
111   fTotNeutralEtAcc = 0;
112   fTotChargedEt  = 0;
113   fTotChargedEtAcc = 0;
114   fMultiplicity = 0;
115   fChargedMultiplicity = 0;
116   fNeutralMultiplicity = 0;
117   
118 }
119 void AliAnalysisHadEt::CreateEtaPtHisto2D(TString name, TString title)
120 {     //creates a 2-d histogram in eta and phi and adds it to the list of histograms to be saved
121   TString *histoname   = new TString();
122   TString *histotitle   = new TString();
123
124   histoname->Append(name);
125   histotitle->Append(title);
126
127   TH2F *histo = new TH2F(histoname->Data(),histotitle->Data(),fgNumOfPtBins, fgPtAxis, fgnumOfEtaBins, fgEtaAxis);
128   histo->SetYTitle("#eta");
129   histo->SetXTitle("p_{T}");
130   histo->SetZTitle("E_{T}");
131   histo->Sumw2();
132   fhistoList->Add(histo);
133   delete histoname;
134   delete histotitle;
135     
136 }
137 void AliAnalysisHadEt::CreateResolutionPtHisto2D(TString name, TString title, TString xtitle, TString ytitle)
138 {     //creates a 2-d histogram in eta and phi and adds it to the list of histograms to be saved
139   TString *histoname   = new TString();
140   TString *histotitle   = new TString();
141
142   histoname->Append(name);
143   histotitle->Append(title);
144
145   //TH2F *histo = new TH2F(histoname->Data(),histotitle->Data(),fgNumOfPtBins, fgPtAxis, fgNumOfResBins, fgResAxis);
146   TH2F *histo = new TH2F(histoname->Data(),histotitle->Data(),fgNumOfPtBins, fgPtAxis, fgNumOfResBins, fgResAxis);
147   histo->SetYTitle(ytitle);
148   histo->SetXTitle(xtitle);
149   histo->SetZTitle("Number of entries");
150   histo->Sumw2();
151   fhistoList->Add(histo);
152   delete histoname;
153   delete histotitle;
154     
155 }
156 void AliAnalysisHadEt::CreatePtHisto1D(TString name, TString title, TString xtitle, TString ytitle)
157 {     //creates a 2-d histogram in eta and phi and adds it to the list of histograms to be saved
158   TString *histoname   = new TString();
159   TString *histotitle   = new TString();
160
161   histoname->Append(name);
162   histotitle->Append(title);
163
164   TH1F *histo = new TH1F(histoname->Data(),histotitle->Data(),fgNumOfPtBins, fgPtAxis);
165   histo->SetYTitle(ytitle);
166   histo->SetXTitle(xtitle);
167   histo->Sumw2();
168   fhistoList->Add(histo);
169   delete histoname;
170   delete histotitle;
171     
172 }
173
174 void AliAnalysisHadEt::CreateHisto1D(TString name, TString title, TString xtitle, TString ytitle,Int_t xbins, Float_t xlow,Float_t xhigh)
175 {     //creates a 1d histogram of the given dimensions and adds it to the list of histograms to be saved
176   TString *histoname   = new TString();
177   TString *histotitle   = new TString();
178   histoname->Append(name);
179   histotitle->Append(title);
180   TH1F *histo = new TH1F(histoname->Data(),histotitle->Data(),xbins,xlow,xhigh);
181   histo->SetYTitle(ytitle);
182   histo->SetXTitle(xtitle);
183   histo->Sumw2();
184   fhistoList->Add(histo);
185   delete histoname;
186   delete histotitle;
187 }
188 void AliAnalysisHadEt::CreateIntHisto1D(TString name, TString title, TString xtitle, TString ytitle,Int_t xbins, Int_t xlow,Int_t xhigh)
189 {     //creates a 1d integer histogram and adds it to the list of histograms to be saved
190   TString *histoname   = new TString();
191   TString *histotitle   = new TString();
192   histoname->Append(name);
193   histotitle->Append(title);
194   TH1I *histo = new TH1I(histoname->Data(),histotitle->Data(),xbins,xlow,xhigh);
195   histo->SetYTitle(ytitle);
196   histo->SetXTitle(xtitle);
197   histo->Sumw2();
198   fhistoList->Add(histo);
199   delete histoname;
200   delete histotitle;
201     
202 }
203 void AliAnalysisHadEt::CreateHisto2D(TString name, TString title, TString xtitle, TString ytitle,Int_t xbins, Float_t xlow,Float_t xhigh,Int_t ybins,Float_t ylow,Float_t yhigh)
204 {     //creates a 2d histogram and adds it to the list of histograms to be saved
205   TString *histoname   = new TString();
206   TString *histotitle   = new TString();
207   histoname->Append(name);
208   histotitle->Append(title);
209   TH2F *histo = new TH2F(histoname->Data(),histotitle->Data(),xbins,xlow,xhigh,ybins,ylow,yhigh);
210   histo->SetYTitle(ytitle);
211   histo->SetXTitle(xtitle);
212   histo->Sumw2();
213   fhistoList->Add(histo);
214   delete histoname;
215   delete histotitle;
216     
217 }
218 void AliAnalysisHadEt::CreateIntHisto2D(TString name, TString title, TString xtitle, TString ytitle,Int_t xbins, Int_t xlow,Int_t xhigh,Int_t ybins,Int_t ylow,Int_t yhigh)
219 {     //creates a 2-d integer histogram and adds it to the list of histograms to be saved
220   TString *histoname   = new TString();
221   TString *histotitle   = new TString();
222   histoname->Append(name);
223   histotitle->Append(title);
224   TH2I *histo = new TH2I(histoname->Data(),histotitle->Data(),xbins,xlow,xhigh,ybins,ylow,yhigh);
225   histo->SetYTitle(ytitle);
226   histo->SetXTitle(xtitle);
227   histo->Sumw2();
228   fhistoList->Add(histo);
229   delete histoname;
230   delete histotitle;
231     
232 }
233
234 void AliAnalysisHadEt::CreateEtaHisto1D(TString name, TString title)
235 {     //creates 1d histogram in eta and adds it to the list of histograms to be saved
236   TString *histoname   = new TString();
237   TString *histotitle   = new TString();
238   histoname->Append(name);
239   histotitle->Append(title);
240   TH1F *histo = new TH1F(histoname->Data(),histotitle->Data(),fgnumOfEtaBins, fgEtaAxis);
241   histo->SetYTitle("E_{T}");
242   histo->SetXTitle("#eta");
243   histo->Sumw2();
244   fhistoList->Add(histo);
245   delete histoname;
246   delete histotitle;
247 }
248 void AliAnalysisHadEt::FillHisto1D(TString histname, Float_t x, Float_t weight)
249 {//fills a 1d histogram with the name histoname with the value x and the weight "weight"
250   TH1F     *histo; 
251   TString  *name   = new TString();
252   name->Append(histname);       
253   histo = (TH1F *)fhistoList->FindObject(name->Data()); 
254   if(histo){
255     histo->Fill((Double_t)x, weight);
256   }
257   else{cerr<<"CorrelationMaker::FillHisto1D: no histogram "<< name->Data()<<endl;}
258   delete name;
259 }
260 void AliAnalysisHadEt::FillHisto2D(TString histname, Float_t x, Float_t y, Float_t weight)
261 {//fills a 2d histogram with the name histoname with the value x and the weight "weight"
262   TH2F     *histo; 
263   TString  *name   = new TString();
264   name->Append(histname);       
265   histo = (TH2F *)fhistoList->FindObject(name->Data()); 
266   if(histo){
267     histo->Fill((Double_t)x,(Double_t)y, weight);
268   }
269   else{cerr<<"CorrelationMaker::FillHisto2D: no histogram "<< name->Data()<<endl;}
270   delete name;
271 }
272
273
274 Float_t AliAnalysisHadEt::Et(TParticle *part, float mass){//function to calculate et in the same way as it would be calculated in a calorimeter
275   if(mass+1000<0.01){//if no mass given return default.  The default argument is -1000
276     if(TMath::Abs(part->GetPDG(0)->PdgCode())==2212 || TMath::Abs(part->GetPDG(0)->PdgCode())==2112){
277       if(part->GetPDG(0)->PdgCode()==-2212 || part->GetPDG(0)->PdgCode()==-2112){//antiproton or antineutron
278         //for antinucleons we specifically want to return the kinetic energy plus twice the rest mass
279         return (part->Energy()+part->GetMass())*TMath::Sin(part->Theta());
280       }
281       if(part->GetPDG(0)->PdgCode()==2212 || part->GetPDG(0)->PdgCode()==2112){//proton or neutron
282         //for nucleons we specifically want to return the kinetic energy only
283         return (part->Energy()-part->GetMass())*TMath::Sin(part->Theta());
284       }
285     }
286     else{//otherwise go to the default
287       return part->Energy()*TMath::Sin(part->Theta());
288     }
289   }
290   else{//otherwise use the mass that was given
291     return (TMath::Sqrt(TMath::Power(part->P(),2.0)+TMath::Power(mass,2.0)))*TMath::Sin(part->Theta());
292   }
293   return 0.0;
294 }
295 Float_t AliAnalysisHadEt::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
296   if(pid==AliAnalysisEtCommon::fgPiPlusCode || pid==AliAnalysisEtCommon::fgPiMinusCode){//Nothing special for pions
297     return TMath::Sqrt(p*p + fgPionMass*fgPionMass) * TMath::Sin(theta);
298   }
299   if(pid==AliAnalysisEtCommon::fgKPlusCode || pid==AliAnalysisEtCommon::fgKMinusCode){//Nothing special for kaons
300     return TMath::Sqrt(p*p + AliAnalysisEtCommon::fgKaonMass*AliAnalysisEtCommon::fgKaonMass) * TMath::Sin(theta);
301   }
302   if(pid==AliAnalysisEtCommon::fgEPlusCode || pid==AliAnalysisEtCommon::fgEMinusCode){//Nothing special for electrons
303     return TMath::Sqrt(p*p + AliAnalysisEtCommon::fgElectronMass*AliAnalysisEtCommon::fgElectronMass) * TMath::Sin(theta);
304   }
305   if(pid==AliAnalysisEtCommon::fgProtonCode || pid==AliAnalysisEtCommon::fgAntiProtonCode){//But for protons we must be careful...
306     if(charge<0.0){//antiprotns: kinetic energy plus twice the rest mass
307       return (TMath::Sqrt(p*p + AliAnalysisEtCommon::fgProtonMass*AliAnalysisEtCommon::fgProtonMass) + AliAnalysisEtCommon::fgProtonMass) * TMath::Sin(theta);
308     }
309     if(charge>0.0){//antiprotns: kinetic energy only
310       return (TMath::Sqrt(p*p + AliAnalysisEtCommon::fgProtonMass*AliAnalysisEtCommon::fgProtonMass) - AliAnalysisEtCommon::fgProtonMass) * TMath::Sin(theta);
311     }
312   }
313   //cerr<<"Uh-oh!  Et not set properly!"<<endl;
314   return 0.0;
315 }
316
317 Float_t AliAnalysisHadEt::TrueP(float pTrec) const {
318   if(pTrec>1.0) return pTrec;
319   return pTrec/(1-599.334*pTrec+7285.15*pTrec*pTrec)+pTrec;
320 }
321
322 Int_t AliAnalysisHadEt::GetCentralityBin(Int_t numberofbins,AliCentrality *centrality){
323   Int_t centralitybin = -1;
324   if(numberofbins<21) centralitybin= centrality->GetCentralityClass10(fCentralityMethod);
325   else{
326     if(numberofbins<41) centralitybin= centrality->GetCentralityClass5(fCentralityMethod);
327     else{
328       Float_t centpercent = centrality->GetCentralityPercentile(fCentralityMethod);
329       centralitybin= centrality->GetCentralityClass5(fCentralityMethod);
330       Float_t centralitybinfloat = (centpercent/2.5);
331       if(centralitybin>=0){
332         centralitybin =(Int_t) (centralitybinfloat);
333         //cout<<" centbin "<<centralitybin<<" centrality "<<centpercent<<" "<<centpercent/2.5 <<endl;
334       }
335     }
336   }
337   //cout<<" centrality bin "<<centralitybin<<endl;
338   return centralitybin;
339 }