1 //_________________________________________________________________________
2 // Utility Class for transverse energy studies
4 // - reconstruction and MC output
7 //*-- Authors: Oystein Djuvsland (Bergen), David Silvermyr (ORNL)
8 //_________________________________________________________________________
9 //Necessary to read config macros
12 #include <TInterpreter.h>
18 #include "AliESDEvent.h"
19 #include "AliMCEvent.h"
20 #include "AliESDtrackCuts.h"
22 #include "AliAnalysisTaskTotEt.h"
23 #include "AliAnalysisEtReconstructedPhos.h"
24 #include "AliAnalysisEtReconstructedEmcal.h"
25 #include "AliAnalysisEtMonteCarloPhos.h"
26 #include "AliAnalysisEtMonteCarloEmcal.h"
32 ClassImp(AliAnalysisTaskTotEt)
34 //________________________________________________________________________
35 AliAnalysisTaskTotEt::AliAnalysisTaskTotEt(const char *name, Bool_t isMc) :
36 AliAnalysisTaskTransverseEnergy(name, isMc)
42 // select if we should use EMCal or PHOS class
43 // PHOS by default, EMCal if name string contains EMC
46 if (t.Contains("EMC")) {
47 if (fMCConfigFile.Length()) {
48 cout<<"Rereading AliAnalysisEtMonteCarloEmcal configuration file..."<<endl;
49 gROOT->LoadMacro(fMCConfigFile);
50 fMCAnalysis = (AliAnalysisEtMonteCarloEmcal *) gInterpreter->ProcessLine("ConfigEtMonteCarlo()");
53 if (fRecoConfigFile.Length()) {
54 cout<<"Rereading AliAnalysisEtReconstructedEmcal configuration file..."<<endl;
55 gROOT->LoadMacro(fRecoConfigFile);
56 fRecAnalysis = (AliAnalysisEtReconstructedEmcal *) gInterpreter->ProcessLine("ConfigEtReconstructed()");
60 if (fMCConfigFile.Length()) {
61 cout<<"Rereading AliAnalysisEtMonteCarloPhos configuration file..."<<endl;
62 gROOT->LoadMacro(fMCConfigFile);
63 fMCAnalysis = (AliAnalysisEtMonteCarloPhos *) gInterpreter->ProcessLine("ConfigEtMonteCarlo()");
66 if (fRecoConfigFile.Length()) {
67 cout<<"Rereading AliAnalysisEtReconstructedPhos configuration file..."<<endl;
68 gROOT->LoadMacro(fRecoConfigFile);
69 fRecAnalysis = (AliAnalysisEtReconstructedPhos *) gInterpreter->ProcessLine("ConfigEtReconstructed()");
72 // Define input and output slots here
73 // Input slot #0 works with a TChain
74 DefineInput(0, TChain::Class());
75 // Output slot #1 writes into a TH1 container
77 DefineOutput(1, TList::Class());
80 AliAnalysisTaskTotEt::~AliAnalysisTaskTotEt() {//Destructor
86 //________________________________________________________________________
87 void AliAnalysisTaskTotEt::UserCreateOutputObjects()
91 fMCAnalysis->CreateHistograms();
92 fRecAnalysis->CreateHistograms();
93 fOutputList = new TList;
94 fOutputList->SetOwner();
95 fRecAnalysis->FillOutputList(fOutputList);
96 fMCAnalysis->FillOutputList(fOutputList);
97 fHistEtRecvsEtMC = new TH2F("fHistEtRecvsEtMC", "Reconstructed E_{t} vs MC E_{t}", 1000, 0.000, 100, 1000, 0.0001, 100);
98 fOutputList->Add(fHistEtRecvsEtMC);
100 Bool_t selectPrimaries=kTRUE;
101 if (fRecAnalysis->DataSet()==2009) {
102 cout<<"Setting track cuts for the 2009 p+p collisions at 900 GeV"<<endl;
103 fEsdtrackCutsITSTPC = AliESDtrackCuts::GetStandardITSTPCTrackCuts2009(selectPrimaries);
104 fEsdtrackCutsITSTPC->SetName("fEsdTrackCuts");
105 fEsdtrackCutsTPC = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
106 fEsdtrackCutsTPC->SetName("fEsdTrackCutsTPCOnly");
107 //ITS stand alone cuts - similar to 2009 cuts but with only ITS hits required
108 fEsdtrackCutsITS = AliESDtrackCuts::GetStandardITSPureSATrackCuts2009(kTRUE,kFALSE);//we do want primaries but we do not want to require PID info
109 fEsdtrackCutsITS->SetName("fEsdTrackCutsITS");
111 if (fRecAnalysis->DataSet()==2010) {
112 //cout<<"Setting track cuts for the 2010 p+p collisions at 7 GeV"<<endl;
113 cout<<"Warning: Have not set 2010 track cuts yet!!"<<endl;
114 fEsdtrackCutsITSTPC = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010(selectPrimaries);
115 fEsdtrackCutsITSTPC->SetName("fEsdTrackCuts");
116 fEsdtrackCutsTPC = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
117 fEsdtrackCutsTPC->SetName("fEsdTrackCutsTPCOnly");
118 //ITS stand alone cuts - similar to 2009 cuts but with only ITS hits required
119 fEsdtrackCutsITS = AliESDtrackCuts::GetStandardITSPureSATrackCuts2009(kTRUE,kFALSE);//we do want primaries but we do not want to require PID info
120 fEsdtrackCutsITS->SetName("fEsdTrackCutsITS");
123 fOutputList->Add(fEsdtrackCutsITSTPC);
124 fOutputList->Add(fEsdtrackCutsTPC);
125 fOutputList->Add(fEsdtrackCutsITS);
126 if (fEsdtrackCutsITSTPC && fEsdtrackCutsTPC) {
127 fRecAnalysis->SetITSTrackCuts( GetITSTrackCuts());
128 fMCAnalysis->SetITSTrackCuts( GetITSTrackCuts());
129 fRecAnalysis->SetTPCITSTrackCuts( GetTPCITSTrackCuts());
130 fMCAnalysis->SetTPCITSTrackCuts( GetTPCITSTrackCuts());
131 fRecAnalysis->SetTPCOnlyTrackCuts( GetTPCOnlyTrackCuts());
132 fMCAnalysis->SetTPCOnlyTrackCuts( GetTPCOnlyTrackCuts());
136 Printf("Error: no track cuts!");
141 //________________________________________________________________________
142 void AliAnalysisTaskTotEt::UserExec(Option_t *)
145 fESDEvent = dynamic_cast<AliESDEvent*>(InputEvent());
148 Printf("ERROR: Could not retrieve event");
152 Int_t res = CheckPhysicsSelection(fESDEvent->GetRunNumber());
154 AliCentrality *cent = GetCentralityObject();
156 if (res == 0 && cent)
158 if (IsPhysicsSelected())
160 fRecAnalysis->SetCentralityObject(cent);
161 fRecAnalysis->AnalyseEvent(fESDEvent);
163 AliMCEvent* mcEvent = MCEvent();
166 fMCAnalysis->AnalyseEvent(mcEvent);
168 fHistEtRecvsEtMC->Fill(fRecAnalysis->GetTotEtAcc(), fMCAnalysis->GetTotEt());
172 PostData(1, fOutputList);
175 //________________________________________________________________________
176 void AliAnalysisTaskTotEt::Terminate(Option_t *)
178 // Draw result to the screen
179 // Called once at the end of the query
181 fOutputList = dynamic_cast<TList*> (GetOutputData(1));
183 printf("ERROR: Output list not available\n");