]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/totEt/AliAnalysisTaskHadEt.cxx
Adding infrastructure for pPb and 8 TeV pp analyses
[u/mrichter/AliRoot.git] / PWGLF / totEt / AliAnalysisTaskHadEt.cxx
1 //_________________________________________________________________________
2 //  Utility Class for transverse energy studies; charged hadrons
3 //  Task for analysis
4 //  - reconstruction and MC output
5 // implementation file
6 //
7 //Created by Christine Nattrass, Rebecca Scott, Irakli Martashvili
8 //University of Tennessee at Knoxville
9 //_________________________________________________________________________
10 //Necessary to read config macros
11 #include <TROOT.h>
12 #include <TSystem.h>
13 #include <TInterpreter.h>
14
15 #include "TChain.h"
16 #include "TList.h"
17 #include "TH2F.h"
18
19 #include "AliESDEvent.h"
20 #include "AliMCEvent.h"
21 #include "AliESDtrackCuts.h"
22
23 #include "AliAnalysisTaskHadEt.h"
24 #include "AliAnalysisHadEtReconstructed.h"
25 #include "AliAnalysisHadEtMonteCarlo.h"
26 #include "AliPWG0Helper.h"
27 #include "AliTriggerAnalysis.h"
28 #include "AliAnalysisManager.h"
29 #include "AliPIDResponse.h"
30 #include "AliTPCPIDResponse.h" 
31 #include "AliInputEventHandler.h"
32
33 #include <iostream>
34 #include "AliLog.h"
35
36 using namespace std;
37
38 ClassImp(AliAnalysisTaskHadEt)
39
40
41
42 //________________________________________________________________________
43   AliAnalysisTaskHadEt::AliAnalysisTaskHadEt(const char *name, Bool_t isMc, TString recoConfigFile, TString mcConfigFile) :
44         AliAnalysisTaskTransverseEnergy(name, isMc)
45         ,fPIDResponse(0)
46         ,fRecAnalysis(0)
47         ,fMCAnalysis(0)
48         ,fIsSim(isMc)
49         ,kIsOfflineV0AND(0)
50         ,kIsOfflineMB(0)
51 {
52     // Constructor
53   //input hander
54   AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
55   if (!man) {
56     AliFatal("Analysis manager needed");
57     return;
58   }
59
60   AliInputEventHandler *inputHandler=dynamic_cast<AliInputEventHandler*>(man->GetInputEventHandler());
61   if (!inputHandler) {
62     AliFatal("Input handler needed");
63     return;
64   }
65
66   //pid response object
67   fPIDResponse=inputHandler->GetPIDResponse();
68   if (!fPIDResponse) AliError("PIDResponse object was not created");
69   else{cout<<"PIDResponse was created!"<<endl;}
70
71
72   fMCConfigFile = mcConfigFile;
73   fRecoConfigFile = recoConfigFile;
74
75   if(fMCAnalysis) delete fMCAnalysis;
76   if(fRecAnalysis) delete fRecAnalysis;
77
78   if (fRecoConfigFile.Length()) {
79     cout<<"Rereading AliAnalysisHadEtReconstructed configuration file "<<fRecoConfigFile<<endl;
80     gROOT->LoadMacro(fRecoConfigFile);
81     fRecAnalysis = (AliAnalysisHadEtReconstructed *) gInterpreter->ProcessLine("ConfigHadEtReconstructed()");
82   }
83
84   if (fMCConfigFile.Length()) {
85     cout<<"Rereading AliAnalysisHadEtMonteCarlo configuration file "<<fMCConfigFile<<endl;
86     gROOT->LoadMacro(fMCConfigFile);
87     fMCAnalysis = (AliAnalysisHadEtMonteCarlo *) gInterpreter->ProcessLine("ConfigHadEtMonteCarlo()");
88     fMCAnalysis->SetHadEtReconstructed(fRecAnalysis);
89   }
90
91     // Define input and output slots here
92     // Input slot #0 works with a TChain
93     DefineInput(0, TChain::Class());
94     // Output slot #1 writes into a TH1 container
95
96     DefineOutput(1, TList::Class());
97 }
98 AliAnalysisTaskHadEt::~AliAnalysisTaskHadEt(){//Destructor
99   delete fRecAnalysis;
100   delete fMCAnalysis;
101   delete fPIDResponse;
102 }
103
104
105 //________________________________________________________________________
106 void AliAnalysisTaskHadEt::UserCreateOutputObjects()
107 {
108     // Create histograms
109
110     // Called once
111
112
113   //input hander
114   AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
115   if (!man) {
116     AliFatal("Analysis manager needed");
117     return;
118   }
119
120   AliInputEventHandler *inputHandler=dynamic_cast<AliInputEventHandler*>(man->GetInputEventHandler());
121   if (!inputHandler) {
122     AliFatal("Input handler needed");
123     return;
124   }
125
126   //pid response object
127   fPIDResponse=inputHandler->GetPIDResponse();
128   if (!fPIDResponse) AliError("PIDResponse object was not created");
129   else{cout<<"PIDResponse was created!"<<endl;}
130
131
132   fOutputList = new TList;
133   fOutputList->SetOwner();
134   fMCAnalysis->SetHistoList(fOutputList);
135   fRecAnalysis->SetHistoList(fOutputList);
136   if(fIsSim) fMCAnalysis->CreateHistograms();
137   fRecAnalysis->CreateHistograms();
138
139
140   if(fRecAnalysis->DataSet() != fMCAnalysis->DataSet()){
141     cout<<"Warning: Reconstruction data set and Monte Carlo data set are not the same!  Setting data set to "<<fRecAnalysis->DataSet()<<endl;
142   }
143
144   Bool_t selectPrimaries=kTRUE;
145   if(fEsdtrackCutsITSTPC) delete fEsdtrackCutsITSTPC;
146   if(fEsdtrackCutsITS) delete fEsdtrackCutsITS;
147   if(fEsdtrackCutsTPC) delete fEsdtrackCutsTPC;
148   //We do not use these because we are using the 2010 900 GeV data
149 //   if(fRecAnalysis->DataSet()==2009){
150 //     cout<<"Setting track cuts for the 2009 p+p collisions at 900 GeV"<<endl;
151 //     fEsdtrackCutsITSTPC = AliESDtrackCuts::GetStandardITSTPCTrackCuts2009(selectPrimaries);
152 //     fEsdtrackCutsITSTPC->SetName("fEsdTrackCuts");
153 //     fEsdtrackCutsTPC = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
154 //     fEsdtrackCutsTPC->SetName("fEsdTrackCutsTPCOnly");
155 //     //ITS stand alone cuts - similar to 2009 cuts but with only ITS hits required
156 //     fEsdtrackCutsITS =  AliESDtrackCuts::GetStandardITSPureSATrackCuts2009(kTRUE,kFALSE);//we do want primaries but we do not want to require PID info
157 //     fEsdtrackCutsITS->SetName("fEsdTrackCutsITS");
158 //   }
159   if(fRecAnalysis->DataSet()==2010 || fRecAnalysis->DataSet()==20111||fRecAnalysis->DataSet()==2009 || fRecAnalysis->DataSet()==2012 || fRecAnalysis->DataSet()==2013){
160     // AliAnalysisTaskSE::      SelectCollisionCandidates(AliVEvent::kINT7 ) ;
161     if(fRecAnalysis->DataSet()==2010)cout<<"Setting track cuts for the 2010 p+p collisions at 7 TeV"<<endl;
162     else{
163       if(fRecAnalysis->DataSet()==2012)cout<<"Setting track cuts for the 2012 p+p collisions at 8 TeV"<<endl;
164       else{
165         if(fRecAnalysis->DataSet()==2013)cout<<"Setting track cuts for the 2013 p+Pb collisions at 5 TeV"<<endl;
166         else{
167           if(fRecAnalysis->DataSet()==2009){cout<<"Setting track cuts for the 2010 p+p collisions at 900 GeV"<<endl;}
168           else{cout<<"Setting track cuts for the 2011 p+p collisions at 2.76 TeV"<<endl;}
169         }
170       }
171     }
172     //cout<<"Warning:  Have not set 2010 track cuts yet!!"<<endl;
173     fEsdtrackCutsITSTPC = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010(selectPrimaries);
174     fEsdtrackCutsITSTPC->SetName("fEsdTrackCuts");
175     fEsdtrackCutsTPC = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
176     fEsdtrackCutsTPC->SetName("fEsdTrackCutsTPCOnly");
177     //ITS stand alone cuts - similar to 2009 cuts but with only ITS hits required
178     fEsdtrackCutsITS =  AliESDtrackCuts::GetStandardITSSATrackCuts2010(kTRUE,kFALSE);//we do want primaries but we do not want to require PID info
179     fEsdtrackCutsITS->SetName("fEsdTrackCutsITS");
180   }
181   if(fRecAnalysis->DataSet()==20100){
182     cout<<"Setting track cuts for the 2010 Pb+Pb collisions at 2.76 TeV"<<endl;
183     //cout<<"Warning:  Have not set 2010 track cuts yet!!"<<endl;
184     fEsdtrackCutsITSTPC = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010(selectPrimaries);
185     fEsdtrackCutsITSTPC->SetName("fEsdTrackCuts");
186     fEsdtrackCutsTPC = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
187     fEsdtrackCutsTPC->SetName("fEsdTrackCutsTPCOnly");
188     //ITS stand alone cuts - similar to 2009 cuts but with only ITS hits required
189     fEsdtrackCutsITS =  AliESDtrackCuts::GetStandardITSSATrackCutsPbPb2010(kTRUE,kFALSE);//we do want primaries but we do not want to require PID info
190     // fEsdtrackCutsITS =  AliESDtrackCuts::GetStandardITSPureSATrackCuts2010(kTRUE,kFALSE);//we do want primaries but we do not want to require PID info
191    fEsdtrackCutsITS->SetName("fEsdTrackCutsITS");
192   }
193
194   fOutputList->Add(fEsdtrackCutsITSTPC);
195   fOutputList->Add(fEsdtrackCutsTPC);
196   fOutputList->Add(fEsdtrackCutsITS);
197   if(fEsdtrackCutsITSTPC && fEsdtrackCutsTPC){
198     fRecAnalysis->SetITSTrackCuts( GetITSTrackCuts());
199     fMCAnalysis->SetITSTrackCuts( GetITSTrackCuts());
200     fRecAnalysis->SetTPCITSTrackCuts( GetTPCITSTrackCuts());
201     fMCAnalysis->SetTPCITSTrackCuts( GetTPCITSTrackCuts());
202     fRecAnalysis->SetTPCOnlyTrackCuts( GetTPCOnlyTrackCuts());
203     fMCAnalysis->SetTPCOnlyTrackCuts( GetTPCOnlyTrackCuts());
204     //add ITS stuff!
205   }
206   else{
207     Printf("Error: no track cuts!");
208   }
209
210
211
212  PostData(1, fOutputList);
213 }
214
215 //________________________________________________________________________
216 void AliAnalysisTaskHadEt::UserExec(Option_t *)
217 { // execute method
218   fESDEvent = dynamic_cast<AliESDEvent*>(InputEvent());
219 if (!fESDEvent) {
220   Printf("ERROR: Could not retrieve event");
221   return;
222  }
223 //cout<<"AliAnalysisTaskHadEt 165"<<endl;
224
225 //Int_t res = CheckPhysicsSelection(fESDEvent->GetRunNumber()); // Check if the physics selection is valid for this run
226
227 //AliCentrality *cent = GetCentralityObject();
228
229 //if(res == 0 && cent){
230 //if(cent){
231  AliTriggerAnalysis *fTriggerAnalysis = new AliTriggerAnalysis();
232
233   kIsOfflineV0AND = fTriggerAnalysis->IsOfflineTriggerFired(fESDEvent, AliTriggerAnalysis::kV0AND);  
234   kIsOfflineMB = fTriggerAnalysis->IsOfflineTriggerFired(fESDEvent, AliTriggerAnalysis::kMB1);  
235   fRecAnalysis->SetIsOfflineV0AND(kIsOfflineV0AND);
236   fMCAnalysis->SetIsOfflineV0AND(kIsOfflineV0AND);
237   fMCAnalysis->SetIsOfflineMB(kIsOfflineMB);
238
239   Int_t eventtype =     AliPWG0Helper::kInvalidProcess;
240   if(fIsSim && fRecAnalysis->DataSet()!=20100) eventtype = (Int_t) AliPWG0Helper::GetEventProcessType(MCEvent()->Header());
241   fRecAnalysis->AnalyseEvent(fESDEvent,eventtype);
242   if(fIsSim){
243     AliMCEvent* mcEvent = MCEvent();
244     if(!mcEvent){  
245       AliFatal("ERROR: MC Event does not exist");
246       return;
247     }
248     if (fESDEvent){
249       ((AliAnalysisHadEtMonteCarlo*)fMCAnalysis)->AnalyseEvent((AliVEvent*)mcEvent,(AliVEvent*)fESDEvent);
250
251     if( fRecAnalysis->DataSet()==20100 || AliPWG0Helper::GetEventProcessType(mcEvent->Header()) == AliPWG0Helper::kND){//either non-diffractive or Pb+Pb
252       if(fMCAnalysis->Full()){
253         fMCAnalysis->FillSimTotEtMinusRecoTotEtFullAcceptanceTPC( fRecAnalysis->GetCorrectedTotEtFullAcceptanceTPC() );
254         fMCAnalysis->FillSimTotEtMinusRecoTotEtFullAcceptanceITS( fRecAnalysis->GetCorrectedTotEtFullAcceptanceITS() );
255         fMCAnalysis->FillSimTotEtMinusRecoTotEtFullAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedTotEtFullAcceptanceTPCNoPID() );
256         fMCAnalysis->FillSimTotEtMinusRecoTotEtFullAcceptanceITSNoPID( fRecAnalysis->GetCorrectedTotEtFullAcceptanceITSNoPID() );
257         fMCAnalysis->FillSimHadEtMinusRecoHadEtFullAcceptanceTPC( fRecAnalysis->GetCorrectedHadEtFullAcceptanceTPC() );
258         fMCAnalysis->FillSimHadEtMinusRecoHadEtFullAcceptanceITS( fRecAnalysis->GetCorrectedHadEtFullAcceptanceITS() );
259         fMCAnalysis->FillSimHadEtMinusRecoHadEtFullAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedHadEtFullAcceptanceTPCNoPID() );
260         fMCAnalysis->FillSimHadEtMinusRecoHadEtFullAcceptanceITSNoPID( fRecAnalysis->GetCorrectedHadEtFullAcceptanceITSNoPID() );
261
262         fMCAnalysis->FillSimTotEtVsRecoTotEtFullAcceptanceTPC( fRecAnalysis->GetCorrectedTotEtFullAcceptanceTPC() );
263         fMCAnalysis->FillSimTotEtVsRecoTotEtFullAcceptanceITS( fRecAnalysis->GetCorrectedTotEtFullAcceptanceITS() );
264         fMCAnalysis->FillSimTotEtVsRecoTotEtFullAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedTotEtFullAcceptanceTPCNoPID() );
265         fMCAnalysis->FillSimTotEtVsRecoTotEtFullAcceptanceITSNoPID( fRecAnalysis->GetCorrectedTotEtFullAcceptanceITSNoPID() );
266         fMCAnalysis->FillSimHadEtVsRecoHadEtFullAcceptanceTPC( fRecAnalysis->GetCorrectedHadEtFullAcceptanceTPC() );
267         fMCAnalysis->FillSimHadEtVsRecoHadEtFullAcceptanceITS( fRecAnalysis->GetCorrectedHadEtFullAcceptanceITS() );
268         fMCAnalysis->FillSimHadEtVsRecoHadEtFullAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedHadEtFullAcceptanceTPCNoPID() );
269         fMCAnalysis->FillSimHadEtVsRecoHadEtFullAcceptanceITSNoPID( fRecAnalysis->GetCorrectedHadEtFullAcceptanceITSNoPID() );
270
271         fMCAnalysis->FillSimPiKPEtVsRecoPiKPEtFullAcceptanceTPC( fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceTPC() );
272         fMCAnalysis->FillSimPiKPEtVsRecoPiKPEtFullAcceptanceITS( fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceITS() );
273         fMCAnalysis->FillSimPiKPEtVsRecoPiKPEtFullAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceTPCNoPID() );
274         fMCAnalysis->FillSimPiKPEtVsRecoPiKPEtFullAcceptanceITSNoPID( fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceITSNoPID() );//Had
275       }
276       if(fMCAnalysis->EMCAL()){
277         fMCAnalysis->FillSimTotEtMinusRecoTotEtEMCALAcceptanceTPC( fRecAnalysis->GetCorrectedTotEtEMCALAcceptanceTPC() );
278         fMCAnalysis->FillSimTotEtMinusRecoTotEtEMCALAcceptanceITS( fRecAnalysis->GetCorrectedTotEtEMCALAcceptanceITS() );
279         fMCAnalysis->FillSimTotEtMinusRecoTotEtEMCALAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedTotEtEMCALAcceptanceTPCNoPID() );
280         fMCAnalysis->FillSimTotEtMinusRecoTotEtEMCALAcceptanceITSNoPID( fRecAnalysis->GetCorrectedTotEtEMCALAcceptanceITSNoPID() );
281         fMCAnalysis->FillSimHadEtMinusRecoHadEtEMCALAcceptanceTPC( fRecAnalysis->GetCorrectedHadEtEMCALAcceptanceTPC() );
282         fMCAnalysis->FillSimHadEtMinusRecoHadEtEMCALAcceptanceITS( fRecAnalysis->GetCorrectedHadEtEMCALAcceptanceITS() );
283         fMCAnalysis->FillSimHadEtMinusRecoHadEtEMCALAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedHadEtEMCALAcceptanceTPCNoPID() );
284         fMCAnalysis->FillSimHadEtMinusRecoHadEtEMCALAcceptanceITSNoPID( fRecAnalysis->GetCorrectedHadEtEMCALAcceptanceITSNoPID() );
285       }
286       if(fMCAnalysis->PHOS()){
287         fMCAnalysis->FillSimTotEtMinusRecoTotEtPHOSAcceptanceTPC( fRecAnalysis->GetCorrectedTotEtPHOSAcceptanceTPC() );
288         fMCAnalysis->FillSimTotEtMinusRecoTotEtPHOSAcceptanceITS( fRecAnalysis->GetCorrectedTotEtPHOSAcceptanceITS() );
289         fMCAnalysis->FillSimTotEtMinusRecoTotEtPHOSAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedTotEtPHOSAcceptanceTPCNoPID() );
290         fMCAnalysis->FillSimTotEtMinusRecoTotEtPHOSAcceptanceITSNoPID( fRecAnalysis->GetCorrectedTotEtPHOSAcceptanceITSNoPID() );
291         fMCAnalysis->FillSimHadEtMinusRecoHadEtPHOSAcceptanceTPC( fRecAnalysis->GetCorrectedHadEtPHOSAcceptanceTPC() );
292         fMCAnalysis->FillSimHadEtMinusRecoHadEtPHOSAcceptanceITS( fRecAnalysis->GetCorrectedHadEtPHOSAcceptanceITS() );
293         fMCAnalysis->FillSimHadEtMinusRecoHadEtPHOSAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedHadEtPHOSAcceptanceTPCNoPID() );
294         fMCAnalysis->FillSimHadEtMinusRecoHadEtPHOSAcceptanceITSNoPID( fRecAnalysis->GetCorrectedHadEtPHOSAcceptanceITSNoPID() );
295       }
296       if(fMCAnalysis->PiKP() && fMCAnalysis->Full()){
297         fMCAnalysis->FillSimPiKPMinusRecoPiKPFullAcceptanceTPC(fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceTPC());
298         fMCAnalysis->FillSimPiKPMinusRecoPiKPFullAcceptanceITS(fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceITS());
299         fMCAnalysis->FillSimPiKPMinusRecoPiKPFullAcceptanceTPCNoPID(fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceTPCNoPID());
300         fMCAnalysis->FillSimPiKPMinusRecoPiKPFullAcceptanceITSNoPID(fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceITSNoPID());
301       }
302     }
303   }
304   }
305   delete fTriggerAnalysis;
306   //}
307 //cout<<"End Event"<<endl<<endl;
308 // Post output data.
309  PostData(1, fOutputList);
310 }
311
312 //________________________________________________________________________
313 void AliAnalysisTaskHadEt::Terminate(Option_t *)
314 {
315     // Draw result to the screen
316     // Called once at the end of the query
317
318     fOutputList = dynamic_cast<TList*> (GetOutputData(1));
319     if (!fOutputList) {
320         printf("ERROR: Output list not available\n");
321         return;
322     }
323 }
324
325