]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG4/totEt/AliAnalysisHadEtReconstructed.cxx
Changing range on PbPb histograms
[u/mrichter/AliRoot.git] / PWG4 / totEt / AliAnalysisHadEtReconstructed.cxx
CommitLineData
cf6522d1 1//_________________________________________________________________________
2// Utility Class for transverse energy studies, charged hadrons
3// Base class for ESD analysis
4// - reconstruction output
5// implementation file
6//
7//Created by Christine Nattrass, Rebecca Scott, Irakli Martashvili
641e1e0c 8//University of Tennessee at Knoxville
cf6522d1 9//_________________________________________________________________________
3ce6b879 10
11#include <TROOT.h>
12#include <TSystem.h>
13#include <TInterpreter.h>
641e1e0c 14#include "AliAnalysisHadEtReconstructed.h"
15#include "AliAnalysisEtCuts.h"
16#include "AliESDtrack.h"
17#include "AliESDCaloCluster.h"
641e1e0c 18#include "AliVEvent.h"
19#include "AliESDEvent.h"
3ce6b879 20#include "AliESDtrackCuts.h"
21#include "AliESDpid.h"
641e1e0c 22#include "AliVParticle.h"
23#include <iostream>
3ce6b879 24#include "AliAnalysisHadEtCorrections.h"
25#include "TFile.h"
26#include "TString.h"
7d2d1773 27#include "AliAnalysisEtCommon.h"
28#include "AliAnalysisHadEt.h"
0f70cf50 29#include "AliCentrality.h"
0f6416f3 30#include "AliLog.h"
641e1e0c 31
16abb579 32using namespace std;
33
34ClassImp(AliAnalysisHadEtReconstructed);
35
36
641e1e0c 37AliAnalysisHadEtReconstructed::AliAnalysisHadEtReconstructed() :
0f70cf50 38 AliAnalysisHadEt()
39 ,fCorrections(0)
40 ,fConfigFile("ConfigHadEtAnalysis.C")
41 ,fCorrTotEtFullAcceptanceTPC(0)
42 ,fCorrTotEtFullAcceptanceITS(0)
43 ,fCorrHadEtFullAcceptanceTPC(0)
44 ,fCorrHadEtFullAcceptanceITS(0)
45 ,fCorrTotEtEMCALAcceptanceTPC(0)
46 ,fCorrTotEtEMCALAcceptanceITS(0)
47 ,fCorrHadEtEMCALAcceptanceTPC(0)
48 ,fCorrHadEtEMCALAcceptanceITS(0)
49 ,fCorrTotEtPHOSAcceptanceTPC(0)
50 ,fCorrTotEtPHOSAcceptanceITS(0)
51 ,fCorrHadEtPHOSAcceptanceTPC(0)
52 ,fCorrHadEtPHOSAcceptanceITS(0)
53 ,fCorrectedHadEtFullAcceptanceTPCNoPID(0)
54 ,fCorrectedHadEtFullAcceptanceITSNoPID(0)
55 ,fCorrectedHadEtEMCALAcceptanceTPCNoPID(0)
56 ,fCorrectedHadEtEMCALAcceptanceITSNoPID(0)
57 ,fCorrectedHadEtPHOSAcceptanceTPCNoPID(0)
58 ,fCorrectedHadEtPHOSAcceptanceITSNoPID(0)
59 ,fCorrectedHadEtFullAcceptanceTPC(0)
60 ,fCorrectedHadEtFullAcceptanceITS(0)
61 ,fCorrectedHadEtFullAcceptanceTPCAssumingPion(0)
62 ,fCorrectedHadEtFullAcceptanceITSAssumingPion(0)
63 ,fCorrectedHadEtFullAcceptanceTPCAssumingProton(0)
64 ,fCorrectedHadEtFullAcceptanceITSAssumingProton(0)
65 ,fCorrectedHadEtFullAcceptanceTPCAssumingKaon(0)
66 ,fCorrectedHadEtFullAcceptanceITSAssumingKaon(0)
67 ,fCorrectedHadEtEMCALAcceptanceTPC(0)
68 ,fCorrectedHadEtEMCALAcceptanceITS(0)
69 ,fCorrectedHadEtPHOSAcceptanceTPC(0)
70 ,fCorrectedHadEtPHOSAcceptanceITS(0)
71 ,fRawEtFullAcceptanceTPC(0)
72 ,fRawEtFullAcceptanceITS(0)
73 ,fRawEtEMCALAcceptanceTPC(0)
74 ,fRawEtEMCALAcceptanceITS(0)
75 ,fRawEtPHOSAcceptanceTPC(0)
76 ,fRawEtPHOSAcceptanceITS(0)
77 ,fRawEtFullAcceptanceTPCNoPID(0)
78 ,fRawEtFullAcceptanceITSNoPID(0)
79 ,fRawEtEMCALAcceptanceTPCNoPID(0)
80 ,fRawEtEMCALAcceptanceITSNoPID(0)
81 ,fRawEtPHOSAcceptanceTPCNoPID(0)
82 ,fRawEtPHOSAcceptanceITSNoPID(0)
641e1e0c 83{
641e1e0c 84}
85
cf6522d1 86AliAnalysisHadEtReconstructed::~AliAnalysisHadEtReconstructed()
641e1e0c 87{
464aa50c 88 delete fCorrections;
cf6522d1 89}
90
91Int_t AliAnalysisHadEtReconstructed::AnalyseEvent(AliVEvent* ev)
92{ // analyse ESD event
0f70cf50 93 ResetEventValues();
6a0df78a 94 if(!ev){
0f6416f3 95 AliFatal("ERROR: Event does not exist");
6a0df78a 96 return 0;
97 }
3ce6b879 98
0f70cf50 99 AliESDEvent *realEvent = dynamic_cast<AliESDEvent*>(ev);
100 if(!realEvent){
101 AliFatal("ERROR: ESD Event does not exist");
102 return 0;
103 }
104 fCentBin= -1;
a180fac9 105 fGoodEvent = kTRUE;//for p+p collisions if we made it this far we have a good event
0f70cf50 106 if(fDataSet==20100){//If this is Pb+Pb
107 AliCentrality *centrality = realEvent->GetCentrality();
a180fac9 108 if(fNCentBins<21) fCentBin= centrality->GetCentralityClass10(fCentralityMethod);
0f70cf50 109 else{ fCentBin= centrality->GetCentralityClass5(fCentralityMethod);}
a180fac9 110 if(fCentBin ==-1) fGoodEvent = kFALSE;//but for Pb+Pb events we don't want to count events where we did not find a centrality
0f70cf50 111 }
112 //for PID
113 AliESDpid *pID = new AliESDpid();
114 pID->MakePID(realEvent);
115 TString *strTPC = new TString("TPC");
116 TString *strITS = new TString("ITS");
117 TString *strTPCITS = new TString("TPCITS");
118 for(Int_t cutset=0;cutset<2;cutset++){
119 bool isTPC = false;
120 TString *cutName = NULL;
121 TObjArray* list = NULL;
122 switch(cutset){
123 case 0:
124 cutName = strTPCITS;
125 list = fEsdtrackCutsITSTPC->GetAcceptedTracks(realEvent);
126 isTPC = true;
127 break;
128 case 1:
129 cutName = strITS;
130 list = fEsdtrackCutsITS->GetAcceptedTracks(realEvent);
131 break;
132 case 2:
133 cutName = strTPC;
134 list = fEsdtrackCutsTPC->GetAcceptedTracks(realEvent);
135 break;
136 default:
137 cerr<<"Error: cannot fill histograms!"<<endl;
138 return -1;
ec956c46 139 }
0f70cf50 140 Int_t nGoodTracks = list->GetEntries();
141 for (Int_t iTrack = 0; iTrack < nGoodTracks; iTrack++)
142 {
641e1e0c 143
641e1e0c 144
0f70cf50 145 AliESDtrack *track = dynamic_cast<AliESDtrack*> (list->At(iTrack));
146 if (!track)
147 {
148 Printf("ERROR: Could not get track %d", iTrack);
149 continue;
150 }
151 else{
152 if(TMath::Abs(track->Eta())>fCorrections->GetEtaCut()) continue;
153 Float_t nSigmaPion,nSigmaProton,nSigmaKaon,nSigmaElectron;
49b25059 154 pID->MakeTPCPID(track);
155 pID->MakeITSPID(track);
0f70cf50 156 if(cutset!=1){
157 nSigmaPion = TMath::Abs(pID->NumberOfSigmasTPC(track,AliPID::kPion));
158 nSigmaProton = TMath::Abs(pID->NumberOfSigmasTPC(track,AliPID::kProton));
159 nSigmaKaon = TMath::Abs(pID->NumberOfSigmasTPC(track,AliPID::kKaon));
160 nSigmaElectron = TMath::Abs(pID->NumberOfSigmasTPC(track,AliPID::kElectron));
161 }
3ce6b879 162 else{
0f70cf50 163 nSigmaPion = TMath::Abs(pID->NumberOfSigmasITS(track,AliPID::kPion));
164 nSigmaProton = TMath::Abs(pID->NumberOfSigmasITS(track,AliPID::kProton));
165 nSigmaKaon = TMath::Abs(pID->NumberOfSigmasITS(track,AliPID::kKaon));
166 nSigmaElectron = TMath::Abs(pID->NumberOfSigmasITS(track,AliPID::kElectron));
167 }
168 // bool isPion = (nSigmaPion<3.0 && nSigmaProton>2.0 && nSigmaKaon>2.0);
169 // bool isElectron = (nSigmaElectron<2.0 && nSigmaPion>4.0 && nSigmaProton>3.0 && nSigmaKaon>3.0);
170 // bool isKaon = (nSigmaPion>3.0 && nSigmaProton>2.0 && nSigmaKaon<2.0);
171 // bool isProton = (nSigmaPion>3.0 && nSigmaProton<2.0 && nSigmaKaon>2.0);
6a0df78a 172 bool isPion = (nSigmaPion<3.0 && nSigmaProton>2.0 && nSigmaKaon>2.0);
173 bool isElectron = (nSigmaElectron<2.0 && nSigmaPion>4.0 && nSigmaProton>3.0 && nSigmaKaon>3.0);
701b4621 174 bool isKaon = (nSigmaPion>3.0 && nSigmaProton>3.0 && nSigmaKaon<3.0 && track->Pt()<0.45);
175 bool isProton = (nSigmaPion>3.0 && nSigmaProton<3.0 && nSigmaKaon>3.0 && track->Pt()<0.9);
641e1e0c 176
0f70cf50 177 bool unidentified = (!isProton && !isKaon && !isElectron && !isPion);
49b25059 178 if(cutset==1){//ITS dE/dx identification requires tighter cuts on the tracks and we don't gain much from that so we won't do it
179 unidentified = true;
180 isPion=false;
181 isElectron=false;
182 isKaon=false;
183 isProton=false;
184 }
0f70cf50 185 Float_t dEdx = track->GetTPCsignal();
186 if(cutset==1) dEdx = track->GetITSsignal();
187 FillHisto2D(Form("dEdxDataAll%s",cutName->Data()),track->P(),dEdx,1.0);
641e1e0c 188
0f70cf50 189 bool inPHOS = IsInPHOS(track);
190 bool inEMCAL = IsInEMCAL(track);
a02dfa56 191
0f70cf50 192 Float_t corrBkgd=0.0;
193 Float_t corrNotID=0.0;
49b25059 194 Float_t corrNoID = fCorrections->GetNotIDCorrectionNoPID(track->Pt());
0f70cf50 195 Float_t corrEff = 0.0;
196 Float_t corrEffNoID = 0.0;
49b25059 197 if(cutset!=1){//TPC
0f70cf50 198 corrBkgd = fCorrections->GetBackgroundCorrectionTPC(track->Pt());
199 corrEffNoID = fCorrections->GetTPCEfficiencyCorrectionHadron(track->Pt(),fCentBin);
200 corrNotID = fCorrections->GetNotIDConstCorrectionTPC();
201 corrNoID = fCorrections->GetNotIDConstCorrectionTPCNoID();
202 }
203 if(cutset==1){//ITS
204 corrBkgd = fCorrections->GetBackgroundCorrectionITS(track->Pt());
205 corrEffNoID = fCorrections->GetITSEfficiencyCorrectionHadron(track->Pt(),fCentBin);
206 corrNotID = fCorrections->GetNotIDConstCorrectionITS();
207 corrNoID = fCorrections->GetNotIDConstCorrectionITSNoID();
208 }
49b25059 209 if(fDataSet==20100){
210 FillHisto2D("fbkgdVsCentralityBin",fCentBin,corrBkgd,1.0);
211 FillHisto2D("fnotIDVsCentralityBin",fCentBin,corrNotID,1.0);
212 FillHisto2D("fpTcutVsCentralityBin",fCentBin,fCorrections->GetpTCutCorrectionTPC(),1.0);
213 if(fCorrHadEtFullAcceptanceTPC>0.0) FillHisto2D("fneutralVsCentralityBin",fCentBin,1.0/fCorrHadEtFullAcceptanceTPC,1.0);
214 if(fCorrections->GetNeutralCorrection()>0.0) FillHisto2D("ConstantCorrectionsVsCentralityBin",fCentBin,1.0/fCorrections->GetNeutralCorrection(),1.0);
215 }
0f70cf50 216 Float_t et = 0.0;
217 Float_t etNoID = Et(track->P(),track->Theta(),fgPiPlusCode,track->Charge());
218 Float_t etpartialcorrected = 0.0;
219 Float_t etpartialcorrectedPion = 0.0;
220 Float_t etpartialcorrectedKaon = 0.0;
221 Float_t etpartialcorrectedProton = 0.0;
222 Float_t etpartialcorrectedNoID = corrNoID*corrBkgd*corrEffNoID*etNoID;
223 FillHisto2D(Form("EtDataRaw%sNoID",cutName->Data()),track->Pt(),track->Eta(),etpartialcorrectedNoID);
3ce6b879 224
0f70cf50 225 if(isPion){
226 FillHisto2D(Form("dEdxDataPion%s",cutName->Data()),track->P(),dEdx,1.0);
227 et = Et(track->P(),track->Theta(),fgPiPlusCode,track->Charge());
49b25059 228 if(cutset==0){corrEff = fCorrections->GetTPCEfficiencyCorrectionPion(track->Pt(),fCentBin);}
0f70cf50 229 etpartialcorrected = et*corrBkgd*corrEff*corrNotID;
49b25059 230 if(corrEff>0.0&&fDataSet==20100)FillHisto2D("feffPionVsCentralityBin",fCentBin,1.0/corrEff,1.0);
0f70cf50 231 if(track->Charge()>0.0){
232 FillHisto2D(Form("EtDataRaw%sPiPlus",cutName->Data()),track->Pt(),track->Eta(),et);
233 FillHisto2D(Form("EtDataCorrected%sPiPlus",cutName->Data()),track->Pt(),track->Eta(),etpartialcorrected);
3ce6b879 234 }
0f70cf50 235 else{
236 FillHisto2D(Form("EtDataRaw%sPiMinus",cutName->Data()),track->Pt(),track->Eta(),et);
237 FillHisto2D(Form("EtDataCorrected%sPiMinus",cutName->Data()),track->Pt(),track->Eta(),etpartialcorrected);
3ce6b879 238 }
0f70cf50 239 }
240 if(isKaon){
241 FillHisto2D(Form("dEdxDataKaon%s",cutName->Data()),track->P(),dEdx,1.0);
242 et = Et(track->P(),track->Theta(),fgKPlusCode,track->Charge());
243 if(cutset==0){corrEff = fCorrections->GetTPCEfficiencyCorrectionKaon(track->Pt(),fCentBin);}
0f70cf50 244 etpartialcorrected = et*corrBkgd*corrEff*corrNotID;
49b25059 245 if(corrEff>0.0&&fDataSet==20100)FillHisto2D("feffKaonVsCentralityBin",fCentBin,1.0/corrEff,1.0);
3ce6b879 246
0f70cf50 247 if(track->Charge()>0.0){
248 FillHisto2D(Form("EtDataRaw%sKPlus",cutName->Data()),track->Pt(),track->Eta(),et);
249 FillHisto2D(Form("EtDataCorrected%sKPlus",cutName->Data()),track->Pt(),track->Eta(),etpartialcorrected);
3ce6b879 250 }
0f70cf50 251 else{
252 FillHisto2D(Form("EtDataRaw%sKMinus",cutName->Data()),track->Pt(),track->Eta(),et);
253 FillHisto2D(Form("EtDataCorrected%sKMinus",cutName->Data()),track->Pt(),track->Eta(),etpartialcorrected);
3ce6b879 254 }
0f70cf50 255 }
256 if(isProton){
257 FillHisto2D(Form("dEdxDataProton%s",cutName->Data()),track->P(),dEdx,1.0);
258 et = Et(track->P(),track->Theta(),fgProtonCode,track->Charge());
259 if(cutset==0){corrEff = fCorrections->GetTPCEfficiencyCorrectionProton(track->Pt(),fCentBin);}
0f70cf50 260 etpartialcorrected = et*corrBkgd*corrEff*corrNotID;
49b25059 261 if(corrEff>0.0&&fDataSet==20100)FillHisto2D("feffProtonVsCentralityBin",fCentBin,1.0/corrEff,1.0);
0f70cf50 262
263 if(track->Charge()>0.0){
264 FillHisto2D(Form("EtDataRaw%sProton",cutName->Data()),track->Pt(),track->Eta(),et);
265 FillHisto2D(Form("EtDataCorrected%sProton",cutName->Data()),track->Pt(),track->Eta(),etpartialcorrected);
3ce6b879 266 }
3d1099f3 267 else{
0f70cf50 268 FillHisto2D(Form("EtDataRaw%sAntiProton",cutName->Data()),track->Pt(),track->Eta(),et);
269 FillHisto2D(Form("EtDataCorrected%sAntiProton",cutName->Data()),track->Pt(),track->Eta(),etpartialcorrected);
0e866ddc 270 }
3ce6b879 271 }
0f70cf50 272 if(isElectron){
273 FillHisto2D(Form("dEdxDataProton%s",cutName->Data()),track->P(),dEdx,1.0);
0f70cf50 274 }
275 if(unidentified){
49b25059 276 if(isPion) cerr<<"I should not be here!! AliAnalysisHadEtReconstructed 273"<<endl;
0f70cf50 277 FillHisto2D(Form("dEdxDataUnidentified%s",cutName->Data()),track->P(),dEdx,1.0);
278 et = Et(track->P(),track->Theta(),fgPiPlusCode,track->Charge());
279 Float_t etProton = Et(track->P(),track->Theta(),fgProtonCode,track->Charge());
280 Float_t etKaon = Et(track->P(),track->Theta(),fgKPlusCode,track->Charge());
49b25059 281 if(corrEff>0.0&&fDataSet==20100)FillHisto2D("feffHadronVsCentralityBin",fCentBin,1.0/corrEff,1.0);
0f70cf50 282 etpartialcorrected = et*corrBkgd*corrEffNoID*corrNotID;
283 etpartialcorrectedPion = et*corrBkgd*corrEffNoID;
284 etpartialcorrectedProton = etProton*corrBkgd*corrEffNoID;
285 etpartialcorrectedKaon = etKaon*corrBkgd*corrEffNoID;
0f70cf50 286 FillHisto2D(Form("EtDataCorrected%sUnidentified",cutName->Data()),track->Pt(),track->Eta(),etpartialcorrected);
287 }
288 else{
289 etpartialcorrectedPion = etpartialcorrected;
290 etpartialcorrectedKaon = etpartialcorrected;
291 etpartialcorrectedProton = etpartialcorrected;
292 }
293 if(!isTPC){
294 etpartialcorrected = etpartialcorrectedNoID;//Not using PID for ITS
295 }
296 AddEt(et,etNoID,etpartialcorrected,etpartialcorrectedPion,etpartialcorrectedProton,etpartialcorrectedKaon,etpartialcorrectedNoID,track->Pt(),isTPC,inPHOS,inEMCAL);
3ce6b879 297 }
0f70cf50 298 }
299 delete list;
300 }
a180fac9 301 if(GetCorrectedHadEtFullAcceptanceTPC()>0.0 && fGoodEvent)FillHisto1D("RecoHadEtFullAcceptanceTPC",GetCorrectedHadEtFullAcceptanceTPC(),1.0);
302 if(GetCorrectedTotEtFullAcceptanceTPC()>0.0 && fGoodEvent)FillHisto1D("RecoTotEtFullAcceptanceTPC",GetCorrectedTotEtFullAcceptanceTPC(),1.0);
303 if(GetCorrectedHadEtFullAcceptanceTPCAssumingPion()>0.0 && fGoodEvent)FillHisto1D("RecoHadEtFullAcceptanceTPCAssumingPion",GetCorrectedHadEtFullAcceptanceTPCAssumingPion(),1.0);
304 if(GetCorrectedHadEtFullAcceptanceTPCAssumingProton()>0.0 && fGoodEvent)FillHisto1D("RecoHadEtFullAcceptanceTPCAssumingProton",GetCorrectedHadEtFullAcceptanceTPCAssumingProton(),1.0);
305 if(GetCorrectedHadEtFullAcceptanceTPCAssumingKaon()>0.0 && fGoodEvent)FillHisto1D("RecoHadEtFullAcceptanceTPCAssumingKaon",GetCorrectedHadEtFullAcceptanceTPCAssumingKaon(),1.0);
306 if(GetCorrectedHadEtEMCALAcceptanceTPC()>0.0 && fGoodEvent)FillHisto1D("RecoHadEtEMCALAcceptanceTPC",GetCorrectedHadEtEMCALAcceptanceTPC(),1.0);
307 if(GetCorrectedTotEtEMCALAcceptanceTPC()>0.0 && fGoodEvent)FillHisto1D("RecoTotEtEMCALAcceptanceTPC",GetCorrectedTotEtEMCALAcceptanceTPC(),1.0);
308 if(GetCorrectedHadEtPHOSAcceptanceTPC()>0.0 && fGoodEvent)FillHisto1D("RecoHadEtPHOSAcceptanceTPC",GetCorrectedHadEtPHOSAcceptanceTPC(),1.0);
309 if(GetCorrectedTotEtPHOSAcceptanceTPC()>0.0 && fGoodEvent)FillHisto1D("RecoTotEtPHOSAcceptanceTPC",GetCorrectedTotEtPHOSAcceptanceTPC(),1.0);
310 if(GetCorrectedHadEtFullAcceptanceTPCNoPID()>0.0 && fGoodEvent)FillHisto1D("RecoHadEtFullAcceptanceTPCNoPID",GetCorrectedHadEtFullAcceptanceTPCNoPID(),1.0);
311 if(GetCorrectedTotEtFullAcceptanceTPCNoPID()>0.0 && fGoodEvent)FillHisto1D("RecoTotEtFullAcceptanceTPCNoPID",GetCorrectedTotEtFullAcceptanceTPCNoPID(),1.0);
312 if(GetCorrectedHadEtEMCALAcceptanceTPCNoPID()>0.0 && fGoodEvent)FillHisto1D("RecoHadEtEMCALAcceptanceTPCNoPID",GetCorrectedHadEtEMCALAcceptanceTPCNoPID(),1.0);
313 if(GetCorrectedTotEtEMCALAcceptanceTPCNoPID()>0.0 && fGoodEvent)FillHisto1D("RecoTotEtEMCALAcceptanceTPCNoPID",GetCorrectedTotEtEMCALAcceptanceTPCNoPID(),1.0);
314 if(GetCorrectedHadEtPHOSAcceptanceTPCNoPID()>0.0 && fGoodEvent)FillHisto1D("RecoHadEtPHOSAcceptanceTPCNoPID",GetCorrectedHadEtPHOSAcceptanceTPCNoPID(),1.0);
315 if(GetCorrectedTotEtPHOSAcceptanceTPCNoPID()>0.0 && fGoodEvent)FillHisto1D("RecoTotEtPHOSAcceptanceTPCNoPID",GetCorrectedTotEtPHOSAcceptanceTPCNoPID(),1.0);
316 if(GetCorrectedHadEtFullAcceptanceITS()>0.0 && fGoodEvent)FillHisto1D("RecoHadEtFullAcceptanceITS",GetCorrectedHadEtFullAcceptanceITS(),1.0);
317 if(GetCorrectedHadEtFullAcceptanceITSAssumingPion()>0.0 && fGoodEvent)FillHisto1D("RecoHadEtFullAcceptanceITSAssumingPion",GetCorrectedHadEtFullAcceptanceITSAssumingPion(),1.0);
318 if(GetCorrectedHadEtFullAcceptanceITSAssumingProton()>0.0 && fGoodEvent)FillHisto1D("RecoHadEtFullAcceptanceITSAssumingProton",GetCorrectedHadEtFullAcceptanceITSAssumingProton(),1.0);
319 if(GetCorrectedHadEtFullAcceptanceITSAssumingKaon()>0.0 && fGoodEvent)FillHisto1D("RecoHadEtFullAcceptanceITSAssumingKaon",GetCorrectedHadEtFullAcceptanceITSAssumingKaon(),1.0);
320 if(GetCorrectedTotEtFullAcceptanceITS()>0.0 && fGoodEvent)FillHisto1D("RecoTotEtFullAcceptanceITS",GetCorrectedTotEtFullAcceptanceITS(),1.0);
321 if(GetCorrectedHadEtEMCALAcceptanceITS()>0.0 && fGoodEvent)FillHisto1D("RecoHadEtEMCALAcceptanceITS",GetCorrectedHadEtEMCALAcceptanceITS(),1.0);
322 if(GetCorrectedTotEtEMCALAcceptanceITS()>0.0 && fGoodEvent)FillHisto1D("RecoTotEtEMCALAcceptanceITS",GetCorrectedTotEtEMCALAcceptanceITS(),1.0);
323 if(GetCorrectedHadEtPHOSAcceptanceITS()>0.0 && fGoodEvent)FillHisto1D("RecoHadEtPHOSAcceptanceITS",GetCorrectedHadEtPHOSAcceptanceITS(),1.0);
324 if(GetCorrectedTotEtPHOSAcceptanceITS()>0.0 && fGoodEvent)FillHisto1D("RecoTotEtPHOSAcceptanceITS",GetCorrectedTotEtPHOSAcceptanceITS(),1.0);
325 if(GetCorrectedHadEtFullAcceptanceITSNoPID()>0.0 && fGoodEvent)FillHisto1D("RecoHadEtFullAcceptanceITSNoPID",GetCorrectedHadEtFullAcceptanceITSNoPID(),1.0);
326 if(GetCorrectedTotEtFullAcceptanceITSNoPID()>0.0 && fGoodEvent)FillHisto1D("RecoTotEtFullAcceptanceITSNoPID",GetCorrectedTotEtFullAcceptanceITSNoPID(),1.0);
327 if(GetCorrectedHadEtEMCALAcceptanceITSNoPID()>0.0 && fGoodEvent)FillHisto1D("RecoHadEtEMCALAcceptanceITSNoPID",GetCorrectedHadEtEMCALAcceptanceITSNoPID(),1.0);
328 if(GetCorrectedTotEtEMCALAcceptanceITSNoPID()>0.0 && fGoodEvent)FillHisto1D("RecoTotEtEMCALAcceptanceITSNoPID",GetCorrectedTotEtEMCALAcceptanceITSNoPID(),1.0);
329 if(GetCorrectedHadEtPHOSAcceptanceITSNoPID()>0.0 && fGoodEvent)FillHisto1D("RecoHadEtPHOSAcceptanceITSNoPID",GetCorrectedHadEtPHOSAcceptanceITSNoPID(),1.0);
330 if(GetCorrectedTotEtPHOSAcceptanceITSNoPID()>0.0 && fGoodEvent)FillHisto1D("RecoTotEtPHOSAcceptanceITSNoPID",GetCorrectedTotEtPHOSAcceptanceITSNoPID(),1.0);
331
332 if(GetRawEtFullAcceptanceTPC()>0.0 && fGoodEvent)FillHisto1D("RecoRawEtFullAcceptanceTPC",GetRawEtFullAcceptanceTPC(),1.0);
333 if(GetRawEtEMCALAcceptanceTPC()>0.0 && fGoodEvent)FillHisto1D("RecoRawEtEMCALAcceptanceTPC",GetRawEtEMCALAcceptanceTPC(),1.0);
334 if(GetRawEtPHOSAcceptanceTPC()>0.0 && fGoodEvent)FillHisto1D("RecoRawEtPHOSAcceptanceTPC",GetRawEtPHOSAcceptanceTPC(),1.0);
335 if(GetRawEtFullAcceptanceTPCNoPID()>0.0 && fGoodEvent)FillHisto1D("RecoRawEtFullAcceptanceTPCNoPID",GetRawEtFullAcceptanceTPCNoPID(),1.0);
336 if(GetRawEtEMCALAcceptanceTPCNoPID()>0.0 && fGoodEvent)FillHisto1D("RecoRawEtEMCALAcceptanceTPCNoPID",GetRawEtEMCALAcceptanceTPCNoPID(),1.0);
337 if(GetRawEtPHOSAcceptanceTPCNoPID()>0.0 && fGoodEvent)FillHisto1D("RecoRawEtPHOSAcceptanceTPCNoPID",GetRawEtPHOSAcceptanceTPCNoPID(),1.0);
338 if(GetRawEtFullAcceptanceITS()>0.0 && fGoodEvent)FillHisto1D("RecoRawEtFullAcceptanceITS",GetRawEtFullAcceptanceITS(),1.0);
339 if(GetRawEtEMCALAcceptanceITS()>0.0 && fGoodEvent)FillHisto1D("RecoRawEtEMCALAcceptanceITS",GetRawEtEMCALAcceptanceITS(),1.0);
340 if(GetRawEtPHOSAcceptanceITS()>0.0 && fGoodEvent)FillHisto1D("RecoRawEtPHOSAcceptanceITS",GetRawEtPHOSAcceptanceITS(),1.0);
341 if(GetRawEtFullAcceptanceITSNoPID()>0.0 && fGoodEvent)FillHisto1D("RecoRawEtFullAcceptanceITSNoPID",GetRawEtFullAcceptanceITSNoPID(),1.0);
342 if(GetRawEtEMCALAcceptanceITSNoPID()>0.0 && fGoodEvent)FillHisto1D("RecoRawEtEMCALAcceptanceITSNoPID",GetRawEtEMCALAcceptanceITSNoPID(),1.0);
343 if(GetRawEtPHOSAcceptanceITSNoPID()>0.0 && fGoodEvent)FillHisto1D("RecoRawEtPHOSAcceptanceITSNoPID",GetRawEtPHOSAcceptanceITSNoPID(),1.0);
0f70cf50 344 if(fCentBin>-1){//if we have Pb+Pb and found a centrality bin
a180fac9 345 if(GetCorrectedHadEtFullAcceptanceTPC()>0.0 && fGoodEvent)FillHisto1D(Form("RecoHadEtFullAcceptanceTPCCB%i",fCentBin),GetCorrectedHadEtFullAcceptanceTPC(),1.0);
346 if(GetCorrectedTotEtFullAcceptanceTPC()>0.0 && fGoodEvent)FillHisto1D(Form("RecoTotEtFullAcceptanceTPCCB%i",fCentBin),GetCorrectedTotEtFullAcceptanceTPC(),1.0);
347 if(GetCorrectedHadEtFullAcceptanceITS()>0.0 && fGoodEvent)FillHisto1D(Form("RecoHadEtFullAcceptanceITSCB%i",fCentBin),GetCorrectedHadEtFullAcceptanceITS(),1.0);
348 if(GetCorrectedTotEtFullAcceptanceITS()>0.0 && fGoodEvent)FillHisto1D(Form("RecoTotEtFullAcceptanceITSCB%i",fCentBin),GetCorrectedTotEtFullAcceptanceITS(),1.0);
349 if(GetRawEtFullAcceptanceTPC()>0.0 && fGoodEvent) FillHisto1D(Form("RecoRawEtFullAcceptanceTPCCB%i",fCentBin),GetRawEtFullAcceptanceTPC(),1.0);
350 if(GetRawEtFullAcceptanceITS()>0.0 && fGoodEvent) FillHisto1D(Form("RecoRawEtFullAcceptanceITSCB%i",fCentBin),GetRawEtFullAcceptanceITS(),1.0);
0f70cf50 351 }
352 delete pID;
353 delete strTPC;
354 delete strITS;
355 delete strTPCITS;
49b25059 356 // cout<<"Reconstructed pi/k/p et "<<GetCorrectedPiKPEtFullAcceptanceTPC()<<endl;
0f70cf50 357 return 1;
641e1e0c 358}
3d1099f3 359void AliAnalysisHadEtReconstructed::AddEt(Float_t rawEt, Float_t rawEtNoPID, Float_t corrEt, Float_t corrEtPion, Float_t corrEtProton, Float_t corrEtKaon, Float_t corrEtNoPID, Float_t pt, Bool_t IsTPC, Bool_t InPHOS, Bool_t InEMCAL) {//Adding Et to each of the variables that tracks et event by event
d6214a64 360 if(pt>=AliAnalysisHadEt::fgPtTPCCutOff && IsTPC){//TPC tracks
361 //adding to the raw Et
d6214a64 362 fRawEtFullAcceptanceTPC += rawEt;
363 if(InPHOS)fRawEtPHOSAcceptanceTPC += rawEt;
364 if(InEMCAL)fRawEtEMCALAcceptanceTPC += rawEt;
365 fRawEtFullAcceptanceTPCNoPID += rawEtNoPID;
366 if(InPHOS)fRawEtPHOSAcceptanceTPCNoPID += rawEtNoPID;
367 if(InEMCAL)fRawEtEMCALAcceptanceTPCNoPID += rawEtNoPID;
368 //adding to the corrected Et
49b25059 369 fCorrectedHadEtFullAcceptanceTPC += corrEt;//the pi/k/p et
3d1099f3 370 fCorrectedHadEtFullAcceptanceTPCAssumingPion += corrEtPion;
371 fCorrectedHadEtFullAcceptanceTPCAssumingProton += corrEtProton;
372 fCorrectedHadEtFullAcceptanceTPCAssumingKaon += corrEtKaon;
d6214a64 373 if(InPHOS)fCorrectedHadEtPHOSAcceptanceTPC += corrEt;
374 if(InEMCAL)fCorrectedHadEtEMCALAcceptanceTPC += corrEt;
375 fCorrectedHadEtFullAcceptanceTPCNoPID += corrEtNoPID;
376 if(InPHOS)fCorrectedHadEtPHOSAcceptanceTPCNoPID += corrEtNoPID;
377 if(InEMCAL)fCorrectedHadEtEMCALAcceptanceTPCNoPID += corrEtNoPID;
378 }
379 if(pt<AliAnalysisHadEt::fgPtTPCCutOff &&pt>=AliAnalysisHadEt::fgPtITSCutOff && !IsTPC){//ITS tracks
380 //adding to the raw Et
381 fRawEtFullAcceptanceITS += rawEt;
382 if(InPHOS)fRawEtPHOSAcceptanceITS += rawEt;
383 if(InEMCAL)fRawEtEMCALAcceptanceITS += rawEt;
384 fRawEtFullAcceptanceITSNoPID += rawEtNoPID;
385 if(InPHOS)fRawEtPHOSAcceptanceITSNoPID += rawEtNoPID;
386 if(InEMCAL)fRawEtEMCALAcceptanceITSNoPID += rawEtNoPID;
387 //adding to the corrected Et
388 fCorrectedHadEtFullAcceptanceITS += corrEt;
3d1099f3 389 fCorrectedHadEtFullAcceptanceITSAssumingPion += corrEtPion;
390 fCorrectedHadEtFullAcceptanceITSAssumingProton += corrEtProton;
391 fCorrectedHadEtFullAcceptanceITSAssumingKaon += corrEtKaon;
d6214a64 392 if(InPHOS)fCorrectedHadEtPHOSAcceptanceITS += corrEt;
393 if(InEMCAL)fCorrectedHadEtEMCALAcceptanceITS += corrEt;
394 fCorrectedHadEtFullAcceptanceITSNoPID += corrEtNoPID;
395 if(InPHOS)fCorrectedHadEtPHOSAcceptanceITSNoPID += corrEtNoPID;
396 if(InEMCAL)fCorrectedHadEtEMCALAcceptanceITSNoPID += corrEtNoPID;
397 }
398}
641e1e0c 399
d6214a64 400Bool_t AliAnalysisHadEtReconstructed::IsInPHOS(AliESDtrack *track){//This function will need to be elaborated on later to include PHOS dead channels
6a0df78a 401 if(!track){
402 cout<<"Error: Track does not exist!!"<<endl;
403 return kFALSE;
404 }
d6214a64 405 return TMath::Abs(track->Eta()) < fCuts->GetGeometryPhosEtaAccCut()//in eta acceptance
406 && track->Phi()*180.0/TMath::Pi() > fCuts->GetGeometryPhosPhiAccMinCut()//greater than the minimum phi
407 && track->Phi()*180.0/TMath::Pi() < fCuts->GetGeometryPhosPhiAccMaxCut();//less than the maximum phi
408}
409Bool_t AliAnalysisHadEtReconstructed::IsInEMCAL(AliESDtrack *track){//This function will need to be elaborated on later to include EMCAL dead channels
6a0df78a 410 if(!track){
411 cout<<"Error: Track does not exist!!"<<endl;
412 return kFALSE;
413 }
d6214a64 414 return TMath::Abs(track->Eta()) < fCuts->GetGeometryEmcalEtaAccCut()//in eta acceptance
415 && track->Phi()*180.0/TMath::Pi() > fCuts->GetGeometryEmcalPhiAccMinCut()//greater than the minimum phi
416 && track->Phi()*180.0/TMath::Pi() < fCuts->GetGeometryEmcalPhiAccMaxCut();//less than the maximum phi
417}
418Bool_t AliAnalysisHadEtReconstructed::CheckGoodVertex(AliVParticle* track)
cf6522d1 419{ // check vertex
641e1e0c 420
0f70cf50 421 Float_t bxy = 999.;
422 Float_t bz = 999.;
423 if(!track){
424 AliError("ERROR: no track");
425 return kFALSE;
426 }
427 AliESDtrack *esdTrack = dynamic_cast<AliESDtrack*>(track);
428 if(!esdTrack){
429 AliError("ERROR: no track");
430 return kFALSE;
431 }
432 esdTrack->GetImpactParametersTPC(bxy,bz);
641e1e0c 433
0f70cf50 434 bool status = (TMath::Abs(track->Xv()) < fCuts->GetReconstructedVertexXCut()) &&
435 (TMath::Abs(track->Yv()) < fCuts->GetReconstructedVertexYCut()) &&
436 (TMath::Abs(track->Zv()) < fCuts->GetReconstructedVertexZCut()) &&
437 (TMath::Abs(bxy) < fCuts->GetReconstructedIPxyCut()) &&
438 (TMath::Abs(bz) < fCuts->GetReconstructedIPzCut());
641e1e0c 439
0f70cf50 440 return status;
641e1e0c 441}
442
443void AliAnalysisHadEtReconstructed::Init()
cf6522d1 444{ // Init
4b40b2b1 445 AliAnalysisHadEt::Init();
964c8159 446 if(fCorrections){
464aa50c 447 fCorrTotEtFullAcceptanceTPC = fCorrections->GetConstantCorrections(kTRUE,fgPtTPCCutOff,"Full");
448 fCorrTotEtFullAcceptanceITS = fCorrections->GetConstantCorrections(kTRUE,fgPtITSCutOff,"Full");
449 fCorrHadEtFullAcceptanceTPC = fCorrections->GetConstantCorrections(kFALSE,fgPtTPCCutOff,"Full");
450 fCorrHadEtFullAcceptanceITS = fCorrections->GetConstantCorrections(kFALSE,fgPtITSCutOff,"Full");
451 fCorrTotEtEMCALAcceptanceTPC = fCorrections->GetConstantCorrections(kTRUE,fgPtTPCCutOff,"EMCAL");
452 fCorrTotEtEMCALAcceptanceITS = fCorrections->GetConstantCorrections(kTRUE,fgPtITSCutOff,"EMCAL");
453 fCorrHadEtEMCALAcceptanceTPC = fCorrections->GetConstantCorrections(kFALSE,fgPtTPCCutOff,"EMCAL");
454 fCorrHadEtEMCALAcceptanceITS = fCorrections->GetConstantCorrections(kFALSE,fgPtITSCutOff,"EMCAL");
455 fCorrTotEtPHOSAcceptanceTPC = fCorrections->GetConstantCorrections(kTRUE,fgPtTPCCutOff,"PHOS");
456 fCorrTotEtPHOSAcceptanceITS = fCorrections->GetConstantCorrections(kTRUE,fgPtITSCutOff,"PHOS");
457 fCorrHadEtPHOSAcceptanceTPC = fCorrections->GetConstantCorrections(kFALSE,fgPtTPCCutOff,"PHOS");
458 fCorrHadEtPHOSAcceptanceITS = fCorrections->GetConstantCorrections(kFALSE,fgPtITSCutOff,"PHOS");
964c8159 459 }
460 else{
461 cout<<"Warning! You have not set corrections. Your code will crash. You have to set the corrections."<<endl;
3ce6b879 462 }
641e1e0c 463}
464
464aa50c 465void AliAnalysisHadEtReconstructed::ResetEventValues(){//resetting event by event et's
d6214a64 466 AliAnalysisHadEt::ResetEventValues();
0f70cf50 467 fCorrectedHadEtFullAcceptanceTPCNoPID=0.0;
468 fCorrectedHadEtFullAcceptanceITSNoPID=0.0;
469 fCorrectedHadEtEMCALAcceptanceTPCNoPID=0.0;
470 fCorrectedHadEtEMCALAcceptanceITSNoPID=0.0;
471 fCorrectedHadEtPHOSAcceptanceTPCNoPID=0.0;
472 fCorrectedHadEtPHOSAcceptanceITSNoPID=0.0;
473 fCorrectedHadEtFullAcceptanceTPC=0.0;
474 fCorrectedHadEtFullAcceptanceITS=0.0;
475 fCorrectedHadEtFullAcceptanceTPCAssumingPion=0.0;
476 fCorrectedHadEtFullAcceptanceITSAssumingPion=0.0;
477 fCorrectedHadEtFullAcceptanceTPCAssumingProton=0.0;
478 fCorrectedHadEtFullAcceptanceITSAssumingProton=0.0;
479 fCorrectedHadEtFullAcceptanceTPCAssumingKaon=0.0;
480 fCorrectedHadEtFullAcceptanceITSAssumingKaon=0.0;
481 fCorrectedHadEtEMCALAcceptanceTPC=0.0;
482 fCorrectedHadEtEMCALAcceptanceITS=0.0;
483 fCorrectedHadEtPHOSAcceptanceTPC=0.0;
484 fCorrectedHadEtPHOSAcceptanceITS=0.0;
485 fRawEtFullAcceptanceTPC=0.0;
486 fRawEtFullAcceptanceITS=0.0;
487 fRawEtEMCALAcceptanceTPC=0.0;
488 fRawEtEMCALAcceptanceITS=0.0;
489 fRawEtPHOSAcceptanceTPC=0.0;
490 fRawEtPHOSAcceptanceITS=0.0;
491 fRawEtFullAcceptanceTPCNoPID=0.0;
492 fRawEtFullAcceptanceITSNoPID=0.0;
493 fRawEtEMCALAcceptanceTPCNoPID=0.0;
494 fRawEtEMCALAcceptanceITSNoPID=0.0;
495 fRawEtPHOSAcceptanceTPCNoPID=0.0;
496 fRawEtPHOSAcceptanceITSNoPID=0.0;
d6214a64 497
0f70cf50 498 if(TMath::Abs(fCorrTotEtFullAcceptanceTPC)<1e-3){
499 if (fConfigFile.Length()) {
500 cout<<"Warning: Rereading fCorrections file..."<<endl;
501 gROOT->LoadMacro(fConfigFile);
502 fCorrections = (AliAnalysisHadEtCorrections *) gInterpreter->ProcessLine("ConfigHadEtAnalysis()");
503 fCorrTotEtFullAcceptanceTPC = fCorrections->GetConstantCorrections(kTRUE,fgPtTPCCutOff,"Full");
504 fCorrTotEtFullAcceptanceITS = fCorrections->GetConstantCorrections(kTRUE,fgPtITSCutOff,"Full");
505 fCorrHadEtFullAcceptanceTPC = fCorrections->GetConstantCorrections(kFALSE,fgPtTPCCutOff,"Full");
506 fCorrHadEtFullAcceptanceITS = fCorrections->GetConstantCorrections(kFALSE,fgPtITSCutOff,"Full");
507 fCorrTotEtEMCALAcceptanceTPC = fCorrections->GetConstantCorrections(kTRUE,fgPtTPCCutOff,"EMCAL");
508 fCorrTotEtEMCALAcceptanceITS = fCorrections->GetConstantCorrections(kTRUE,fgPtITSCutOff,"EMCAL");
509 fCorrHadEtEMCALAcceptanceTPC = fCorrections->GetConstantCorrections(kFALSE,fgPtTPCCutOff,"EMCAL");
510 fCorrHadEtEMCALAcceptanceITS = fCorrections->GetConstantCorrections(kFALSE,fgPtITSCutOff,"EMCAL");
511 fCorrTotEtPHOSAcceptanceTPC = fCorrections->GetConstantCorrections(kTRUE,fgPtTPCCutOff,"PHOS");
512 fCorrTotEtPHOSAcceptanceITS = fCorrections->GetConstantCorrections(kTRUE,fgPtITSCutOff,"PHOS");
513 fCorrHadEtPHOSAcceptanceTPC = fCorrections->GetConstantCorrections(kFALSE,fgPtTPCCutOff,"PHOS");
514 fCorrHadEtPHOSAcceptanceITS = fCorrections->GetConstantCorrections(kFALSE,fgPtITSCutOff,"PHOS");
515 }
516 else{cerr<<"Uh-oh! Unable to open configuration file!"<<endl;}
517 }
d6214a64 518}
464aa50c 519void AliAnalysisHadEtReconstructed::CreateHistograms(){//Creating histograms and adding them to the output TList
3ce6b879 520
6a0df78a 521 //TString *strTPC = new TString("TPC");
3ce6b879 522 TString *strITS = new TString("ITS");
523 TString *strTPCITS = new TString("TPCITS");
524 for(Int_t i=0;i<2;i++){
ea331c5d 525 TString *cutName = NULL;
3ce6b879 526 Float_t maxPtdEdx = 10;
527 Float_t mindEdx = 35;
528 Float_t maxdEdx = 150.0;
529 switch(i){
530 case 0:
f43fc416 531 cutName = strTPCITS;
3ce6b879 532 break;
533 case 1:
f43fc416 534 cutName = strITS;
3ce6b879 535 maxPtdEdx = 5;
536 maxdEdx = 500.0;
537 break;
6a0df78a 538 //not deleting this completely since we might want to add it back in later
0f70cf50 539 // case 2:
540 // cutName = strTPC;
541 // maxPtdEdx = 5;
542 // maxdEdx = 500.0;
543 // break;
3ce6b879 544 default:
545 cerr<<"Error: cannot make histograms!"<<endl;
546 return;
547 }
548
549 CreateEtaPtHisto2D(Form("EtDataRaw%sPiPlus",cutName->Data()),"Raw reconstructed E_{T} from identified #pi^{+}");
550 CreateEtaPtHisto2D(Form("EtDataRaw%sPiMinus",cutName->Data()),"Raw reconstructed E_{T} from identified #pi^{-}");
551 CreateEtaPtHisto2D(Form("EtDataRaw%sKPlus",cutName->Data()),"Raw reconstructed E_{T} from identified K^{+}");
0f70cf50 552 // CreateEtaPtHisto2D(Form("EtDataRaw%sEMinus",cutName->Data()),"Raw reconstructed E_{T} from identified e^{-}");
553 // CreateEtaPtHisto2D(Form("EtDataRaw%sEPlus",cutName->Data()),"Raw reconstructed E_{T} from identified e^{+}");
3ce6b879 554 CreateEtaPtHisto2D(Form("EtDataRaw%sKMinus",cutName->Data()),"Raw reconstructed E_{T} from identified K^{-}");
555 CreateEtaPtHisto2D(Form("EtDataRaw%sProton",cutName->Data()),"Raw reconstructed E_{T} from identified p");
556 CreateEtaPtHisto2D(Form("EtDataRaw%sAntiProton",cutName->Data()),"Raw reconstructed E_{T} from identified #bar{p}");
557 CreateEtaPtHisto2D(Form("EtDataRaw%sUnidentified",cutName->Data()),"Raw reconstructed E_{T} from unidentified particles using real mass");
558 CreateEtaPtHisto2D(Form("EtDataRaw%sNoID",cutName->Data()),"Raw reconstructed E_{T} from unidentified particles using real mass");
559
560 CreateEtaPtHisto2D(Form("EtDataCorrected%sPiPlus",cutName->Data()),"Corrected reconstructed E_{T} from identified #pi^{+}");
561 CreateEtaPtHisto2D(Form("EtDataCorrected%sPiMinus",cutName->Data()),"Corrected reconstructed E_{T} from identified #pi^{-}");
562 CreateEtaPtHisto2D(Form("EtDataCorrected%sKPlus",cutName->Data()),"Corrected reconstructed E_{T} from identified K^{+}");
0f70cf50 563 // CreateEtaPtHisto2D(Form("EtDataCorrected%sEMinus",cutName->Data()),"Corrected reconstructed E_{T} from identified e^{-}");
564 // CreateEtaPtHisto2D(Form("EtDataCorrected%sEPlus",cutName->Data()),"Corrected reconstructed E_{T} from identified e^{+}");
3ce6b879 565 CreateEtaPtHisto2D(Form("EtDataCorrected%sKMinus",cutName->Data()),"Corrected reconstructed E_{T} from identified K^{-}");
566 CreateEtaPtHisto2D(Form("EtDataCorrected%sProton",cutName->Data()),"Corrected reconstructed E_{T} from identified p");
567 CreateEtaPtHisto2D(Form("EtDataCorrected%sAntiProton",cutName->Data()),"Corrected reconstructed E_{T} from identified #bar{p}");
568 CreateEtaPtHisto2D(Form("EtDataCorrected%sUnidentified",cutName->Data()),"Corrected reconstructed E_{T} from unidentified particles using real mass");
569 CreateEtaPtHisto2D(Form("EtDataCorrected%sNoID",cutName->Data()),"Corrected reconstructed E_{T} from unidentified particles using real mass");
570
571
572 CreateEtaPtHisto2D(Form("EtNData%sPiPlus",cutName->Data()),"Number of reconstructed #pi^{+}");
573 CreateEtaPtHisto2D(Form("EtNData%sPiMinus",cutName->Data()),"Number of reconstructed #pi^{-}");
574 CreateEtaPtHisto2D(Form("EtNData%sKPlus",cutName->Data()),"Number of reconstructed K^{+}");
575 CreateEtaPtHisto2D(Form("EtNData%sKMinus",cutName->Data()),"Number of reconstructed K^{-}");
576 CreateEtaPtHisto2D(Form("EtNData%sProton",cutName->Data()),"Number of reconstructed p");
577 CreateEtaPtHisto2D(Form("EtNData%sAntiProton",cutName->Data()),"Number of reconstructed #bar{p}");
578 CreateEtaPtHisto2D(Form("EtNData%sUnidentified",cutName->Data()),"Number of Reconstructed unidentified particles");
579
580 CreateHisto2D(Form("dEdxDataAll%s",cutName->Data()),"dE/dx for all particles","momentum (GeV/c)","dE/dx",400,0.0,maxPtdEdx,200,mindEdx,maxdEdx);
581 CreateHisto2D(Form("dEdxDataPion%s",cutName->Data()),"dE/dx for #pi^{#pm}","momentum (GeV/c)","dE/dx",400,0.0,maxPtdEdx,200,mindEdx,maxdEdx);
582 CreateHisto2D(Form("dEdxDataKaon%s",cutName->Data()),"dE/dx for K^{#pm}","momentum (GeV/c)","dE/dx",400,0.0,maxPtdEdx,200,mindEdx,maxdEdx);
583 CreateHisto2D(Form("dEdxDataProton%s",cutName->Data()),"dE/dx for p(#bar{p})","momentum (GeV/c)","dE/dx",400,0.0,maxPtdEdx,200,mindEdx,maxdEdx);
584 CreateHisto2D(Form("dEdxDataElectron%s",cutName->Data()),"dE/dx for e^{#pm}","momentum (GeV/c)","dE/dx",400,0.0,maxPtdEdx,200,mindEdx,maxdEdx);
585 CreateHisto2D(Form("dEdxDataUnidentified%s",cutName->Data()),"dE/dx for unidentified particles","momentum (GeV/c)","dE/dx",400,0.0,maxPtdEdx,200,mindEdx,maxdEdx);
586 }
4e0c0fe1 587
588 Float_t minEt = 0.0;
589 Float_t maxEt = 100.0;
8becf95f 590 if(fDataSet==20100) maxEt=3500.0;
4e0c0fe1 591 Int_t nbinsEt = 200;
592 char histoname[200];
593 char histotitle[200];
594 char xtitle[50];
595 TString *ytitle = new TString("Number of events");
464aa50c 596 TString *sTPC = new TString("TPC");
597 TString *sITS = new TString("ITS");
598 TString *sTPCpt = new TString("0.15");
599 TString *sITSpt = new TString("0.10");
600 TString *sPID = new TString("");
601 TString *sNoPID = new TString("NoPID");
602 TString *sNoPIDString = new TString(", No PID");
603 TString *sHadEt = new TString("HadEt");
604 TString *sRawEt = new TString("RawEt");
605 TString *sTotEt = new TString("TotEt");
606 TString *sTotEtString = new TString("total E_{T}");
607 TString *sHadEtString = new TString("hadronic E_{T}");
608 TString *sRawEtString = new TString("raw E_{T}");
609 TString *sFull = new TString("Full");
610 TString *sEMCAL = new TString("EMCAL");
611 TString *sPHOS = new TString("PHOS");
4e0c0fe1 612
613 for(int tpc = 0;tpc<2;tpc++){
614 for(int hadet = 0;hadet<3;hadet++){
615 for(int type = 0;type<3;type++){
616 for(int pid = 0;pid<2;pid++){
ea331c5d 617 TString *detector = NULL;
618 TString *partid = NULL;
464aa50c 619 TString *et = sHadEt;
ea331c5d 620 TString *acceptance = NULL;
621 TString *ptstring = NULL;
622 TString *partidstring = NULL;
464aa50c 623 TString *etstring = sHadEtString;
624 if(tpc==1) {detector = sTPC; ptstring = sTPCpt;}
625 else{detector = sITS; ptstring = sITSpt;}
626 if(pid==1){partid = sPID; partidstring = sPID;}
627 else{partid = sNoPID; partidstring = sNoPIDString;}
628 if(hadet==1) {et = sHadEt; etstring = sHadEtString;}
629 if(hadet==0){et = sTotEt; etstring = sTotEtString;}
630 if(hadet==2){et = sRawEt; etstring = sRawEtString;}
4e0c0fe1 631 switch(type){
632 case 0:
464aa50c 633 acceptance = sFull;
4e0c0fe1 634 break;
635 case 1:
464aa50c 636 acceptance = sEMCAL;
4e0c0fe1 637 break;
638 case 2:
464aa50c 639 acceptance = sPHOS;
4e0c0fe1 640 break;
641 default:
464aa50c 642 acceptance = sFull;
4e0c0fe1 643 }
6a0df78a 644 snprintf(histoname,200,"Reco%s%sAcceptance%s%s",et->Data(),acceptance->Data(),detector->Data(),partid->Data());
645 snprintf(histotitle,200,"Reconstructed %s with %s acceptance for p_{T}>%s GeV/c%s",etstring->Data(),acceptance->Data(),ptstring->Data(),partidstring->Data());
646 snprintf(xtitle,50,"Reconstructed %s",etstring->Data());
0e866ddc 647 CreateHisto1D(histoname,histotitle,xtitle,ytitle->Data(),nbinsEt*2,minEt,maxEt);
0f70cf50 648 if(fDataSet==20100 && type ==0 &&pid==1){//If this is Pb+Pb and full acceptance with pid
649 Int_t width = 5;
a180fac9 650 if(fNCentBins<21) width = 10;
0f70cf50 651 for(Int_t i=0;i<fNCentBins;i++){
652 snprintf(histoname,200,"Reco%s%sAcceptance%s%sCB%i",et->Data(),acceptance->Data(),detector->Data(),partid->Data(),i);
653 snprintf(histotitle,200,"Reconstructed %s with %s acceptance for p_{T}>%s GeV/c%s for centrality %i-%i",etstring->Data(),acceptance->Data(),ptstring->Data(),partidstring->Data(),i*width,(i+1)*width);
654 snprintf(xtitle,50,"Reconstructed %s",etstring->Data());
655 CreateHisto1D(histoname,histotitle,xtitle,ytitle->Data(),nbinsEt*2,minEt,maxEt);
656 }
657 }
4e0c0fe1 658 }
659 }
660 }
661 }
3d1099f3 662 CreateHisto1D("RecoHadEtFullAcceptanceTPCAssumingPion","Reconstructing E_{T}^{had} with full acceptance for p_{T}>0.15 GeV/c assuming pions","Reconstructed E_{T}^{had}","dN_{eve}/dE_{T}^{had}",nbinsEt*2,minEt,maxEt);
663 CreateHisto1D("RecoHadEtFullAcceptanceTPCAssumingProton","Reconstructing E_{T}^{had} with full acceptance for p_{T}>0.15 GeV/c assuming protons","Reconstructed E_{T}^{had}","dN_{eve}/dE_{T}^{had}",nbinsEt*2,minEt,maxEt);
664 CreateHisto1D("RecoHadEtFullAcceptanceTPCAssumingKaon","Reconstructing E_{T}^{had} with full acceptance for p_{T}>0.15 GeV/c assuming kaons","Reconstructed E_{T}^{had}","dN_{eve}/dE_{T}^{had}",nbinsEt*2,minEt,maxEt);
665 CreateHisto1D("RecoHadEtFullAcceptanceITSAssumingPion","Reconstructing E_{T}^{had} with full acceptance for p_{T}>0.10 GeV/c assuming pions","Reconstructed E_{T}^{had}","dN_{eve}/dE_{T}^{had}",nbinsEt*2,minEt,maxEt);
666 CreateHisto1D("RecoHadEtFullAcceptanceITSAssumingProton","Reconstructing E_{T}^{had} with full acceptance for p_{T}>0.10 GeV/c assuming protons","Reconstructed E_{T}^{had}","dN_{eve}/dE_{T}^{had}",nbinsEt*2,minEt,maxEt);
667 CreateHisto1D("RecoHadEtFullAcceptanceITSAssumingKaon","Reconstructing E_{T}^{had} with full acceptance for p_{T}>0.10 GeV/c assuming kaons","Reconstructed E_{T}^{had}","dN_{eve}/dE_{T}^{had}",nbinsEt*2,minEt,maxEt);
0e866ddc 668
a180fac9 669 //Cross checks that corrections are applied correctly
670 if(fDataSet==20100){
671 CreateHisto2D("fbkgdVsCentralityBin","f_{bkgd} vs centrality bin","centrality bin","f_{bkgd}",21,-1.5,19.5,200,0.7,1.05);//
672 CreateHisto2D("feffPionVsCentralityBin","Pion efficiency vs centrality bin","centrality bin","pion efficiency",21,-1.5,19.5,200,0,1.2);//
673 CreateHisto2D("feffHadronVsCentralityBin","Hadron efficiency vs centrality bin","centrality bin","hadron efficiency",21,-1.5,19.5,200,0,1.2);//
674 CreateHisto2D("feffKaonVsCentralityBin","Kaon efficiency vs centrality bin","centrality bin","kaon efficiency",21,-1.5,19.5,200,0,1.2);//
675 CreateHisto2D("feffProtonVsCentralityBin","Proton efficiency vs centrality bin","centrality bin","proton efficiency",21,-1.5,19.5,200,0,1.2);//
676 CreateHisto2D("fnotIDVsCentralityBin","f_{notID} vs centrality bin","centrality bin","f_{notID}",21,-1.5,19.5,50,0.95,1.05);//
677 CreateHisto2D("fpTcutVsCentralityBin","f_{pTcut} vs centrality bin","centrality bin","f_{pTcut}",21,-1.5,19.5,50,0.95,1.05);
678 CreateHisto2D("fneutralVsCentralityBin","f_{neutral} vs centrality bin","centrality bin","f_{neutral}",21,-1.5,19.5,50,0.5,1.00);
679 CreateHisto2D("ConstantCorrectionsVsCentralityBin","constant corrections vs centrality bin","centrality bin","constant corrections",21,-1.5,19.5,50,0.5,1.00);
680 }
681
464aa50c 682 delete sTPC;
683 delete sITS;
684 delete sTPCpt;
685 delete sITSpt;
686 delete sPID;
687 delete sNoPID;
688 delete sNoPIDString;
689 delete sHadEt;
690 delete sTotEt;
691 delete sTotEtString;
692 delete sHadEtString;
693 delete sFull;
694 delete sEMCAL;
695 delete sPHOS;
4e0c0fe1 696
3ce6b879 697}