1 //_________________________________________________________________________
2 // Utility Class for transverse energy studies; charged hadrons
4 // - reconstruction and MC output
7 //Created by Christine Nattrass, Rebecca Scott, Irakli Martashvili
8 //University of Tennessee at Knoxville
9 //_________________________________________________________________________
10 //Necessary to read config macros
13 #include <TInterpreter.h>
19 #include "AliESDEvent.h"
20 #include "AliMCEvent.h"
21 #include "AliESDtrackCuts.h"
23 #include "AliAnalysisTaskHadEt.h"
24 #include "AliAnalysisHadEtReconstructed.h"
25 #include "AliAnalysisHadEtMonteCarlo.h"
26 #include "AliPWG0Helper.h"
27 #include "AliTriggerAnalysis.h"
34 ClassImp(AliAnalysisTaskHadEt)
38 //________________________________________________________________________
39 AliAnalysisTaskHadEt::AliAnalysisTaskHadEt(const char *name, Bool_t isMc, TString recoConfigFile, TString mcConfigFile) :
40 AliAnalysisTaskTransverseEnergy(name, isMc)
47 fMCConfigFile = mcConfigFile;
48 fRecoConfigFile = recoConfigFile;
50 if(fMCAnalysis) delete fMCAnalysis;
51 if(fRecAnalysis) delete fRecAnalysis;
53 if (fRecoConfigFile.Length()) {
54 cout<<"Rereading AliAnalysisHadEtReconstructed configuration file "<<fRecoConfigFile<<endl;
55 gROOT->LoadMacro(fRecoConfigFile);
56 fRecAnalysis = (AliAnalysisHadEtReconstructed *) gInterpreter->ProcessLine("ConfigHadEtReconstructed()");
59 if (fMCConfigFile.Length()) {
60 cout<<"Rereading AliAnalysisHadEtMonteCarlo configuration file "<<fMCConfigFile<<endl;
61 gROOT->LoadMacro(fMCConfigFile);
62 fMCAnalysis = (AliAnalysisHadEtMonteCarlo *) gInterpreter->ProcessLine("ConfigHadEtMonteCarlo()");
63 fMCAnalysis->SetHadEtReconstructed(fRecAnalysis);
66 // Define input and output slots here
67 // Input slot #0 works with a TChain
68 DefineInput(0, TChain::Class());
69 // Output slot #1 writes into a TH1 container
71 DefineOutput(1, TList::Class());
73 AliAnalysisTaskHadEt::~AliAnalysisTaskHadEt(){//Destructor
79 //________________________________________________________________________
80 void AliAnalysisTaskHadEt::UserCreateOutputObjects()
88 fOutputList = new TList;
89 fOutputList->SetOwner();
90 fMCAnalysis->SetHistoList(fOutputList);
91 fRecAnalysis->SetHistoList(fOutputList);
92 if(fIsSim) fMCAnalysis->CreateHistograms();
93 fRecAnalysis->CreateHistograms();
96 if(fRecAnalysis->DataSet() != fMCAnalysis->DataSet()){
97 cout<<"Warning: Reconstruction data set and Monte Carlo data set are not the same! Setting data set to "<<fRecAnalysis->DataSet()<<endl;
100 Bool_t selectPrimaries=kTRUE;
101 if(fEsdtrackCutsITSTPC) delete fEsdtrackCutsITSTPC;
102 if(fEsdtrackCutsITS) delete fEsdtrackCutsITS;
103 if(fEsdtrackCutsTPC) delete fEsdtrackCutsTPC;
104 //We do not use these because we are using the 2010 900 GeV data
105 // if(fRecAnalysis->DataSet()==2009){
106 // cout<<"Setting track cuts for the 2009 p+p collisions at 900 GeV"<<endl;
107 // fEsdtrackCutsITSTPC = AliESDtrackCuts::GetStandardITSTPCTrackCuts2009(selectPrimaries);
108 // fEsdtrackCutsITSTPC->SetName("fEsdTrackCuts");
109 // fEsdtrackCutsTPC = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
110 // fEsdtrackCutsTPC->SetName("fEsdTrackCutsTPCOnly");
111 // //ITS stand alone cuts - similar to 2009 cuts but with only ITS hits required
112 // fEsdtrackCutsITS = AliESDtrackCuts::GetStandardITSPureSATrackCuts2009(kTRUE,kFALSE);//we do want primaries but we do not want to require PID info
113 // fEsdtrackCutsITS->SetName("fEsdTrackCutsITS");
115 if(fRecAnalysis->DataSet()==2010 || fRecAnalysis->DataSet()==20111||fRecAnalysis->DataSet()==2009){
116 // AliAnalysisTaskSE:: SelectCollisionCandidates(AliVEvent::kINT7 ) ;
117 if(fRecAnalysis->DataSet()==2010)cout<<"Setting track cuts for the 2010 p+p collisions at 7 TeV"<<endl;
119 if(fRecAnalysis->DataSet()==2009){cout<<"Setting track cuts for the 2010 p+p collisions at 900 GeV"<<endl;}
120 else{cout<<"Setting track cuts for the 2011 p+p collisions at 2.76 TeV"<<endl;}
122 //cout<<"Warning: Have not set 2010 track cuts yet!!"<<endl;
123 fEsdtrackCutsITSTPC = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010(selectPrimaries);
124 fEsdtrackCutsITSTPC->SetName("fEsdTrackCuts");
125 fEsdtrackCutsTPC = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
126 fEsdtrackCutsTPC->SetName("fEsdTrackCutsTPCOnly");
127 //ITS stand alone cuts - similar to 2009 cuts but with only ITS hits required
128 fEsdtrackCutsITS = AliESDtrackCuts::GetStandardITSSATrackCuts2010(kTRUE,kFALSE);//we do want primaries but we do not want to require PID info
129 fEsdtrackCutsITS->SetName("fEsdTrackCutsITS");
131 if(fRecAnalysis->DataSet()==20100){
132 cout<<"Setting track cuts for the 2010 Pb+Pb collisions at 2.76 TeV"<<endl;
133 //cout<<"Warning: Have not set 2010 track cuts yet!!"<<endl;
134 fEsdtrackCutsITSTPC = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010(selectPrimaries);
135 fEsdtrackCutsITSTPC->SetName("fEsdTrackCuts");
136 fEsdtrackCutsTPC = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
137 fEsdtrackCutsTPC->SetName("fEsdTrackCutsTPCOnly");
138 //ITS stand alone cuts - similar to 2009 cuts but with only ITS hits required
139 fEsdtrackCutsITS = AliESDtrackCuts::GetStandardITSSATrackCutsPbPb2010(kTRUE,kFALSE);//we do want primaries but we do not want to require PID info
140 // fEsdtrackCutsITS = AliESDtrackCuts::GetStandardITSPureSATrackCuts2010(kTRUE,kFALSE);//we do want primaries but we do not want to require PID info
141 fEsdtrackCutsITS->SetName("fEsdTrackCutsITS");
144 fOutputList->Add(fEsdtrackCutsITSTPC);
145 fOutputList->Add(fEsdtrackCutsTPC);
146 fOutputList->Add(fEsdtrackCutsITS);
147 if(fEsdtrackCutsITSTPC && fEsdtrackCutsTPC){
148 fRecAnalysis->SetITSTrackCuts( GetITSTrackCuts());
149 fMCAnalysis->SetITSTrackCuts( GetITSTrackCuts());
150 fRecAnalysis->SetTPCITSTrackCuts( GetTPCITSTrackCuts());
151 fMCAnalysis->SetTPCITSTrackCuts( GetTPCITSTrackCuts());
152 fRecAnalysis->SetTPCOnlyTrackCuts( GetTPCOnlyTrackCuts());
153 fMCAnalysis->SetTPCOnlyTrackCuts( GetTPCOnlyTrackCuts());
157 Printf("Error: no track cuts!");
160 PostData(1, fOutputList);
163 //________________________________________________________________________
164 void AliAnalysisTaskHadEt::UserExec(Option_t *)
166 fESDEvent = dynamic_cast<AliESDEvent*>(InputEvent());
168 Printf("ERROR: Could not retrieve event");
171 //cout<<"AliAnalysisTaskHadEt 165"<<endl;
173 //Int_t res = CheckPhysicsSelection(fESDEvent->GetRunNumber()); // Check if the physics selection is valid for this run
175 //AliCentrality *cent = GetCentralityObject();
177 //if(res == 0 && cent){
179 AliTriggerAnalysis *fTriggerAnalysis = new AliTriggerAnalysis();
181 kIsOfflineV0AND = fTriggerAnalysis->IsOfflineTriggerFired(fESDEvent, AliTriggerAnalysis::kV0AND);
182 fRecAnalysis->SetIsOfflineV0AND(kIsOfflineV0AND);
183 fMCAnalysis->SetIsOfflineV0AND(kIsOfflineV0AND);
185 Int_t eventtype = AliPWG0Helper::kInvalidProcess;
186 if(fIsSim && fRecAnalysis->DataSet()!=20100) eventtype = (Int_t) AliPWG0Helper::GetEventProcessType(MCEvent()->Header());
187 fRecAnalysis->AnalyseEvent(fESDEvent,eventtype);
189 AliMCEvent* mcEvent = MCEvent();
191 AliFatal("ERROR: MC Event does not exist");
195 ((AliAnalysisHadEtMonteCarlo*)fMCAnalysis)->AnalyseEvent((AliVEvent*)mcEvent,(AliVEvent*)fESDEvent);
197 if( fRecAnalysis->DataSet()==20100 || AliPWG0Helper::GetEventProcessType(mcEvent->Header()) == AliPWG0Helper::kND){//either non-diffractive or Pb+Pb
198 if(fMCAnalysis->Full()){
199 fMCAnalysis->FillSimTotEtMinusRecoTotEtFullAcceptanceTPC( fRecAnalysis->GetCorrectedTotEtFullAcceptanceTPC() );
200 fMCAnalysis->FillSimTotEtMinusRecoTotEtFullAcceptanceITS( fRecAnalysis->GetCorrectedTotEtFullAcceptanceITS() );
201 fMCAnalysis->FillSimTotEtMinusRecoTotEtFullAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedTotEtFullAcceptanceTPCNoPID() );
202 fMCAnalysis->FillSimTotEtMinusRecoTotEtFullAcceptanceITSNoPID( fRecAnalysis->GetCorrectedTotEtFullAcceptanceITSNoPID() );
203 fMCAnalysis->FillSimHadEtMinusRecoHadEtFullAcceptanceTPC( fRecAnalysis->GetCorrectedHadEtFullAcceptanceTPC() );
204 fMCAnalysis->FillSimHadEtMinusRecoHadEtFullAcceptanceITS( fRecAnalysis->GetCorrectedHadEtFullAcceptanceITS() );
205 fMCAnalysis->FillSimHadEtMinusRecoHadEtFullAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedHadEtFullAcceptanceTPCNoPID() );
206 fMCAnalysis->FillSimHadEtMinusRecoHadEtFullAcceptanceITSNoPID( fRecAnalysis->GetCorrectedHadEtFullAcceptanceITSNoPID() );
208 fMCAnalysis->FillSimTotEtVsRecoTotEtFullAcceptanceTPC( fRecAnalysis->GetCorrectedTotEtFullAcceptanceTPC() );
209 fMCAnalysis->FillSimTotEtVsRecoTotEtFullAcceptanceITS( fRecAnalysis->GetCorrectedTotEtFullAcceptanceITS() );
210 fMCAnalysis->FillSimTotEtVsRecoTotEtFullAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedTotEtFullAcceptanceTPCNoPID() );
211 fMCAnalysis->FillSimTotEtVsRecoTotEtFullAcceptanceITSNoPID( fRecAnalysis->GetCorrectedTotEtFullAcceptanceITSNoPID() );
212 fMCAnalysis->FillSimHadEtVsRecoHadEtFullAcceptanceTPC( fRecAnalysis->GetCorrectedHadEtFullAcceptanceTPC() );
213 fMCAnalysis->FillSimHadEtVsRecoHadEtFullAcceptanceITS( fRecAnalysis->GetCorrectedHadEtFullAcceptanceITS() );
214 fMCAnalysis->FillSimHadEtVsRecoHadEtFullAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedHadEtFullAcceptanceTPCNoPID() );
215 fMCAnalysis->FillSimHadEtVsRecoHadEtFullAcceptanceITSNoPID( fRecAnalysis->GetCorrectedHadEtFullAcceptanceITSNoPID() );
217 fMCAnalysis->FillSimPiKPEtVsRecoPiKPEtFullAcceptanceTPC( fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceTPC() );
218 fMCAnalysis->FillSimPiKPEtVsRecoPiKPEtFullAcceptanceITS( fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceITS() );
219 fMCAnalysis->FillSimPiKPEtVsRecoPiKPEtFullAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceTPCNoPID() );
220 fMCAnalysis->FillSimPiKPEtVsRecoPiKPEtFullAcceptanceITSNoPID( fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceITSNoPID() );//Had
222 if(fMCAnalysis->EMCAL()){
223 fMCAnalysis->FillSimTotEtMinusRecoTotEtEMCALAcceptanceTPC( fRecAnalysis->GetCorrectedTotEtEMCALAcceptanceTPC() );
224 fMCAnalysis->FillSimTotEtMinusRecoTotEtEMCALAcceptanceITS( fRecAnalysis->GetCorrectedTotEtEMCALAcceptanceITS() );
225 fMCAnalysis->FillSimTotEtMinusRecoTotEtEMCALAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedTotEtEMCALAcceptanceTPCNoPID() );
226 fMCAnalysis->FillSimTotEtMinusRecoTotEtEMCALAcceptanceITSNoPID( fRecAnalysis->GetCorrectedTotEtEMCALAcceptanceITSNoPID() );
227 fMCAnalysis->FillSimHadEtMinusRecoHadEtEMCALAcceptanceTPC( fRecAnalysis->GetCorrectedHadEtEMCALAcceptanceTPC() );
228 fMCAnalysis->FillSimHadEtMinusRecoHadEtEMCALAcceptanceITS( fRecAnalysis->GetCorrectedHadEtEMCALAcceptanceITS() );
229 fMCAnalysis->FillSimHadEtMinusRecoHadEtEMCALAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedHadEtEMCALAcceptanceTPCNoPID() );
230 fMCAnalysis->FillSimHadEtMinusRecoHadEtEMCALAcceptanceITSNoPID( fRecAnalysis->GetCorrectedHadEtEMCALAcceptanceITSNoPID() );
232 if(fMCAnalysis->PHOS()){
233 fMCAnalysis->FillSimTotEtMinusRecoTotEtPHOSAcceptanceTPC( fRecAnalysis->GetCorrectedTotEtPHOSAcceptanceTPC() );
234 fMCAnalysis->FillSimTotEtMinusRecoTotEtPHOSAcceptanceITS( fRecAnalysis->GetCorrectedTotEtPHOSAcceptanceITS() );
235 fMCAnalysis->FillSimTotEtMinusRecoTotEtPHOSAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedTotEtPHOSAcceptanceTPCNoPID() );
236 fMCAnalysis->FillSimTotEtMinusRecoTotEtPHOSAcceptanceITSNoPID( fRecAnalysis->GetCorrectedTotEtPHOSAcceptanceITSNoPID() );
237 fMCAnalysis->FillSimHadEtMinusRecoHadEtPHOSAcceptanceTPC( fRecAnalysis->GetCorrectedHadEtPHOSAcceptanceTPC() );
238 fMCAnalysis->FillSimHadEtMinusRecoHadEtPHOSAcceptanceITS( fRecAnalysis->GetCorrectedHadEtPHOSAcceptanceITS() );
239 fMCAnalysis->FillSimHadEtMinusRecoHadEtPHOSAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedHadEtPHOSAcceptanceTPCNoPID() );
240 fMCAnalysis->FillSimHadEtMinusRecoHadEtPHOSAcceptanceITSNoPID( fRecAnalysis->GetCorrectedHadEtPHOSAcceptanceITSNoPID() );
242 if(fMCAnalysis->PiKP() && fMCAnalysis->Full()){
243 fMCAnalysis->FillSimPiKPMinusRecoPiKPFullAcceptanceTPC(fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceTPC());
244 fMCAnalysis->FillSimPiKPMinusRecoPiKPFullAcceptanceITS(fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceITS());
245 fMCAnalysis->FillSimPiKPMinusRecoPiKPFullAcceptanceTPCNoPID(fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceTPCNoPID());
246 fMCAnalysis->FillSimPiKPMinusRecoPiKPFullAcceptanceITSNoPID(fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceITSNoPID());
251 delete fTriggerAnalysis;
253 //cout<<"End Event"<<endl<<endl;
255 PostData(1, fOutputList);
258 //________________________________________________________________________
259 void AliAnalysisTaskHadEt::Terminate(Option_t *)
261 // Draw result to the screen
262 // Called once at the end of the query
264 fOutputList = dynamic_cast<TList*> (GetOutputData(1));
266 printf("ERROR: Output list not available\n");