]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/totEt/AliAnalysisTaskHadEt.cxx
faa99b83235c16702d61a5dbe789550b2a9e77ed
[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
28 #include <iostream>
29 #include "AliLog.h"
30
31 using namespace std;
32
33 ClassImp(AliAnalysisTaskHadEt)
34
35
36
37 //________________________________________________________________________
38   AliAnalysisTaskHadEt::AliAnalysisTaskHadEt(const char *name, Bool_t isMc, TString recoConfigFile, TString mcConfigFile) :
39         AliAnalysisTaskTransverseEnergy(name, isMc)
40         ,fRecAnalysis(0)
41         ,fMCAnalysis(0)
42         ,fIsSim(isMc)
43 {
44     // Constructor
45   fMCConfigFile = mcConfigFile;
46   fRecoConfigFile = recoConfigFile;
47
48   if(fMCAnalysis) delete fMCAnalysis;
49   if(fRecAnalysis) delete fRecAnalysis;
50
51   if (fRecoConfigFile.Length()) {
52     cout<<"Rereading AliAnalysisHadEtReconstructed configuration file "<<fRecoConfigFile<<endl;
53     gROOT->LoadMacro(fRecoConfigFile);
54     fRecAnalysis = (AliAnalysisHadEtReconstructed *) gInterpreter->ProcessLine("ConfigHadEtReconstructed()");
55   }
56
57   if (fMCConfigFile.Length()) {
58     cout<<"Rereading AliAnalysisHadEtMonteCarlo configuration file "<<fMCConfigFile<<endl;
59     gROOT->LoadMacro(fMCConfigFile);
60     fMCAnalysis = (AliAnalysisHadEtMonteCarlo *) gInterpreter->ProcessLine("ConfigHadEtMonteCarlo()");
61     fMCAnalysis->SetHadEtReconstructed(fRecAnalysis);
62   }
63
64     // Define input and output slots here
65     // Input slot #0 works with a TChain
66     DefineInput(0, TChain::Class());
67     // Output slot #1 writes into a TH1 container
68
69     DefineOutput(1, TList::Class());
70 }
71 AliAnalysisTaskHadEt::~AliAnalysisTaskHadEt(){//Destructor
72   delete fRecAnalysis;
73   delete fMCAnalysis;
74 }
75
76
77 //________________________________________________________________________
78 void AliAnalysisTaskHadEt::UserCreateOutputObjects()
79 {
80     // Create histograms
81
82     // Called once
83
84
85
86   fOutputList = new TList;
87   fOutputList->SetOwner();
88   fMCAnalysis->SetHistoList(fOutputList);
89   fRecAnalysis->SetHistoList(fOutputList);
90   if(fIsSim) fMCAnalysis->CreateHistograms();
91   fRecAnalysis->CreateHistograms();
92
93
94   if(fRecAnalysis->DataSet() != fMCAnalysis->DataSet()){
95     cout<<"Warning: Reconstruction data set and Monte Carlo data set are not the same!  Setting data set to "<<fRecAnalysis->DataSet()<<endl;
96   }
97
98   Bool_t selectPrimaries=kTRUE;
99   if(fEsdtrackCutsITSTPC) delete fEsdtrackCutsITSTPC;
100   if(fEsdtrackCutsITS) delete fEsdtrackCutsITS;
101   if(fEsdtrackCutsTPC) delete fEsdtrackCutsTPC;
102   //We do not use these because we are using the 2010 900 GeV data
103 //   if(fRecAnalysis->DataSet()==2009){
104 //     cout<<"Setting track cuts for the 2009 p+p collisions at 900 GeV"<<endl;
105 //     fEsdtrackCutsITSTPC = AliESDtrackCuts::GetStandardITSTPCTrackCuts2009(selectPrimaries);
106 //     fEsdtrackCutsITSTPC->SetName("fEsdTrackCuts");
107 //     fEsdtrackCutsTPC = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
108 //     fEsdtrackCutsTPC->SetName("fEsdTrackCutsTPCOnly");
109 //     //ITS stand alone cuts - similar to 2009 cuts but with only ITS hits required
110 //     fEsdtrackCutsITS =  AliESDtrackCuts::GetStandardITSPureSATrackCuts2009(kTRUE,kFALSE);//we do want primaries but we do not want to require PID info
111 //     fEsdtrackCutsITS->SetName("fEsdTrackCutsITS");
112 //   }
113   if(fRecAnalysis->DataSet()==2010 || fRecAnalysis->DataSet()==20111||fRecAnalysis->DataSet()==2009){
114     AliAnalysisTaskSE:: SelectCollisionCandidates(AliVEvent::kINT7 ) ;
115     if(fRecAnalysis->DataSet()==2010)cout<<"Setting track cuts for the 2010 p+p collisions at 7 TeV"<<endl;
116     else{
117       if(fRecAnalysis->DataSet()==2009){cout<<"Setting track cuts for the 2010 p+p collisions at 900 GeV"<<endl;}
118       else{cout<<"Setting track cuts for the 2011 p+p collisions at 2.76 TeV"<<endl;}
119     }
120     //cout<<"Warning:  Have not set 2010 track cuts yet!!"<<endl;
121     fEsdtrackCutsITSTPC = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010(selectPrimaries);
122     fEsdtrackCutsITSTPC->SetName("fEsdTrackCuts");
123     fEsdtrackCutsTPC = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
124     fEsdtrackCutsTPC->SetName("fEsdTrackCutsTPCOnly");
125     //ITS stand alone cuts - similar to 2009 cuts but with only ITS hits required
126     fEsdtrackCutsITS =  AliESDtrackCuts::GetStandardITSPureSATrackCuts2010(kTRUE,kFALSE);//we do want primaries but we do not want to require PID info
127     fEsdtrackCutsITS->SetName("fEsdTrackCutsITS");
128   }
129   if(fRecAnalysis->DataSet()==20100){
130     cout<<"Setting track cuts for the 2010 Pb+Pb collisions at 2.76 TeV"<<endl;
131     //cout<<"Warning:  Have not set 2010 track cuts yet!!"<<endl;
132     fEsdtrackCutsITSTPC = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010(selectPrimaries);
133     fEsdtrackCutsITSTPC->SetName("fEsdTrackCuts");
134     fEsdtrackCutsTPC = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
135     fEsdtrackCutsTPC->SetName("fEsdTrackCutsTPCOnly");
136     //ITS stand alone cuts - similar to 2009 cuts but with only ITS hits required
137     fEsdtrackCutsITS =  AliESDtrackCuts::GetStandardITSSATrackCutsPbPb2010(kTRUE,kFALSE);//we do want primaries but we do not want to require PID info
138     // fEsdtrackCutsITS =  AliESDtrackCuts::GetStandardITSPureSATrackCuts2010(kTRUE,kFALSE);//we do want primaries but we do not want to require PID info
139    fEsdtrackCutsITS->SetName("fEsdTrackCutsITS");
140   }
141
142   fOutputList->Add(fEsdtrackCutsITSTPC);
143   fOutputList->Add(fEsdtrackCutsTPC);
144   fOutputList->Add(fEsdtrackCutsITS);
145   if(fEsdtrackCutsITSTPC && fEsdtrackCutsTPC){
146     fRecAnalysis->SetITSTrackCuts( GetITSTrackCuts());
147     fMCAnalysis->SetITSTrackCuts( GetITSTrackCuts());
148     fRecAnalysis->SetTPCITSTrackCuts( GetTPCITSTrackCuts());
149     fMCAnalysis->SetTPCITSTrackCuts( GetTPCITSTrackCuts());
150     fRecAnalysis->SetTPCOnlyTrackCuts( GetTPCOnlyTrackCuts());
151     fMCAnalysis->SetTPCOnlyTrackCuts( GetTPCOnlyTrackCuts());
152     //add ITS stuff!
153   }
154   else{
155     Printf("Error: no track cuts!");
156   }
157
158  PostData(1, fOutputList);
159 }
160
161 //________________________________________________________________________
162 void AliAnalysisTaskHadEt::UserExec(Option_t *)
163 { // execute method
164   fESDEvent = dynamic_cast<AliESDEvent*>(InputEvent());
165 if (!fESDEvent) {
166   Printf("ERROR: Could not retrieve event");
167   return;
168  }
169 //cout<<"AliAnalysisTaskHadEt 165"<<endl;
170
171 Int_t res = CheckPhysicsSelection(fESDEvent->GetRunNumber()); // Check if the physics selection is valid for this run
172
173 AliCentrality *cent = GetCentralityObject();
174
175 //if(res == 0 && cent){
176 //if(cent){
177   
178   Int_t eventtype =     AliPWG0Helper::kInvalidProcess;
179   if(!fIsSim){
180     fRecAnalysis->AnalyseEvent(fESDEvent,eventtype);
181   }
182   else{
183     AliMCEvent* mcEvent = MCEvent();
184     if(!mcEvent){  
185       AliFatal("ERROR: MC Event does not exist");
186       return;
187     }
188     eventtype = (Int_t) AliPWG0Helper::GetEventProcessType(mcEvent->Header());
189     if (fESDEvent){
190       ((AliAnalysisHadEtMonteCarlo*)fMCAnalysis)->AnalyseEvent((AliVEvent*)mcEvent,(AliVEvent*)fESDEvent);
191
192     if(AliPWG0Helper::GetEventProcessType(mcEvent->Header()) == AliPWG0Helper::kND || fRecAnalysis->DataSet()==20100){//either non-diffractive or Pb+Pb
193       if(fMCAnalysis->Full()){
194         fMCAnalysis->FillSimTotEtMinusRecoTotEtFullAcceptanceTPC( fRecAnalysis->GetCorrectedTotEtFullAcceptanceTPC() );
195         fMCAnalysis->FillSimTotEtMinusRecoTotEtFullAcceptanceITS( fRecAnalysis->GetCorrectedTotEtFullAcceptanceITS() );
196         fMCAnalysis->FillSimTotEtMinusRecoTotEtFullAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedTotEtFullAcceptanceTPCNoPID() );
197         fMCAnalysis->FillSimTotEtMinusRecoTotEtFullAcceptanceITSNoPID( fRecAnalysis->GetCorrectedTotEtFullAcceptanceITSNoPID() );
198         fMCAnalysis->FillSimHadEtMinusRecoHadEtFullAcceptanceTPC( fRecAnalysis->GetCorrectedHadEtFullAcceptanceTPC() );
199         fMCAnalysis->FillSimHadEtMinusRecoHadEtFullAcceptanceITS( fRecAnalysis->GetCorrectedHadEtFullAcceptanceITS() );
200         fMCAnalysis->FillSimHadEtMinusRecoHadEtFullAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedHadEtFullAcceptanceTPCNoPID() );
201         fMCAnalysis->FillSimHadEtMinusRecoHadEtFullAcceptanceITSNoPID( fRecAnalysis->GetCorrectedHadEtFullAcceptanceITSNoPID() );
202
203         fMCAnalysis->FillSimTotEtVsRecoTotEtFullAcceptanceTPC( fRecAnalysis->GetCorrectedTotEtFullAcceptanceTPC() );
204         fMCAnalysis->FillSimTotEtVsRecoTotEtFullAcceptanceITS( fRecAnalysis->GetCorrectedTotEtFullAcceptanceITS() );
205         fMCAnalysis->FillSimTotEtVsRecoTotEtFullAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedTotEtFullAcceptanceTPCNoPID() );
206         fMCAnalysis->FillSimTotEtVsRecoTotEtFullAcceptanceITSNoPID( fRecAnalysis->GetCorrectedTotEtFullAcceptanceITSNoPID() );
207         fMCAnalysis->FillSimHadEtVsRecoHadEtFullAcceptanceTPC( fRecAnalysis->GetCorrectedHadEtFullAcceptanceTPC() );
208         fMCAnalysis->FillSimHadEtVsRecoHadEtFullAcceptanceITS( fRecAnalysis->GetCorrectedHadEtFullAcceptanceITS() );
209         fMCAnalysis->FillSimHadEtVsRecoHadEtFullAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedHadEtFullAcceptanceTPCNoPID() );
210         fMCAnalysis->FillSimHadEtVsRecoHadEtFullAcceptanceITSNoPID( fRecAnalysis->GetCorrectedHadEtFullAcceptanceITSNoPID() );
211
212         fMCAnalysis->FillSimPiKPEtVsRecoPiKPEtFullAcceptanceTPC( fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceTPC() );
213         fMCAnalysis->FillSimPiKPEtVsRecoPiKPEtFullAcceptanceITS( fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceITS() );
214         fMCAnalysis->FillSimPiKPEtVsRecoPiKPEtFullAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceTPCNoPID() );
215         fMCAnalysis->FillSimPiKPEtVsRecoPiKPEtFullAcceptanceITSNoPID( fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceITSNoPID() );//Had
216       }
217       if(fMCAnalysis->EMCAL()){
218         fMCAnalysis->FillSimTotEtMinusRecoTotEtEMCALAcceptanceTPC( fRecAnalysis->GetCorrectedTotEtEMCALAcceptanceTPC() );
219         fMCAnalysis->FillSimTotEtMinusRecoTotEtEMCALAcceptanceITS( fRecAnalysis->GetCorrectedTotEtEMCALAcceptanceITS() );
220         fMCAnalysis->FillSimTotEtMinusRecoTotEtEMCALAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedTotEtEMCALAcceptanceTPCNoPID() );
221         fMCAnalysis->FillSimTotEtMinusRecoTotEtEMCALAcceptanceITSNoPID( fRecAnalysis->GetCorrectedTotEtEMCALAcceptanceITSNoPID() );
222         fMCAnalysis->FillSimHadEtMinusRecoHadEtEMCALAcceptanceTPC( fRecAnalysis->GetCorrectedHadEtEMCALAcceptanceTPC() );
223         fMCAnalysis->FillSimHadEtMinusRecoHadEtEMCALAcceptanceITS( fRecAnalysis->GetCorrectedHadEtEMCALAcceptanceITS() );
224         fMCAnalysis->FillSimHadEtMinusRecoHadEtEMCALAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedHadEtEMCALAcceptanceTPCNoPID() );
225         fMCAnalysis->FillSimHadEtMinusRecoHadEtEMCALAcceptanceITSNoPID( fRecAnalysis->GetCorrectedHadEtEMCALAcceptanceITSNoPID() );
226       }
227       if(fMCAnalysis->PHOS()){
228         fMCAnalysis->FillSimTotEtMinusRecoTotEtPHOSAcceptanceTPC( fRecAnalysis->GetCorrectedTotEtPHOSAcceptanceTPC() );
229         fMCAnalysis->FillSimTotEtMinusRecoTotEtPHOSAcceptanceITS( fRecAnalysis->GetCorrectedTotEtPHOSAcceptanceITS() );
230         fMCAnalysis->FillSimTotEtMinusRecoTotEtPHOSAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedTotEtPHOSAcceptanceTPCNoPID() );
231         fMCAnalysis->FillSimTotEtMinusRecoTotEtPHOSAcceptanceITSNoPID( fRecAnalysis->GetCorrectedTotEtPHOSAcceptanceITSNoPID() );
232         fMCAnalysis->FillSimHadEtMinusRecoHadEtPHOSAcceptanceTPC( fRecAnalysis->GetCorrectedHadEtPHOSAcceptanceTPC() );
233         fMCAnalysis->FillSimHadEtMinusRecoHadEtPHOSAcceptanceITS( fRecAnalysis->GetCorrectedHadEtPHOSAcceptanceITS() );
234         fMCAnalysis->FillSimHadEtMinusRecoHadEtPHOSAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedHadEtPHOSAcceptanceTPCNoPID() );
235         fMCAnalysis->FillSimHadEtMinusRecoHadEtPHOSAcceptanceITSNoPID( fRecAnalysis->GetCorrectedHadEtPHOSAcceptanceITSNoPID() );
236       }
237       if(fMCAnalysis->PiKP() && fMCAnalysis->Full()){
238         fMCAnalysis->FillSimPiKPMinusRecoPiKPFullAcceptanceTPC(fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceTPC());
239         fMCAnalysis->FillSimPiKPMinusRecoPiKPFullAcceptanceITS(fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceITS());
240         fMCAnalysis->FillSimPiKPMinusRecoPiKPFullAcceptanceTPCNoPID(fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceTPCNoPID());
241         fMCAnalysis->FillSimPiKPMinusRecoPiKPFullAcceptanceITSNoPID(fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceITSNoPID());
242       }
243     }
244   }
245   }
246   //}
247 //cout<<"End Event"<<endl<<endl;
248 // Post output data.
249  PostData(1, fOutputList);
250 }
251
252 //________________________________________________________________________
253 void AliAnalysisTaskHadEt::Terminate(Option_t *)
254 {
255     // Draw result to the screen
256     // Called once at the end of the query
257
258     fOutputList = dynamic_cast<TList*> (GetOutputData(1));
259     if (!fOutputList) {
260         printf("ERROR: Output list not available\n");
261         return;
262     }
263 }
264
265