]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/totEt/AliAnalysisTaskTotEt.cxx
AliESDCentrality replaced by AliCentrality
[u/mrichter/AliRoot.git] / PWG4 / totEt / AliAnalysisTaskTotEt.cxx
1 //_________________________________________________________________________
2 //  Utility Class for transverse energy studies
3 //  Task for analysis
4 //  - reconstruction and MC output
5 // implementation file
6 //
7 //*-- Authors: Oystein Djuvsland (Bergen), David Silvermyr (ORNL)
8 //_________________________________________________________________________
9 //Necessary to read config macros
10 #include <TROOT.h>
11 #include <TSystem.h>
12 #include <TInterpreter.h>
13
14 #include "TChain.h"
15 #include "TList.h"
16 #include "TH2F.h"
17
18 #include "AliESDEvent.h"
19 #include "AliMCEvent.h"
20 #include "AliESDtrackCuts.h"
21
22 #include "AliAnalysisTaskTotEt.h"
23 #include "AliAnalysisEtReconstructedPhos.h"
24 #include "AliAnalysisEtReconstructedEmcal.h"
25 #include "AliAnalysisEtMonteCarloPhos.h"
26 #include "AliAnalysisEtMonteCarloEmcal.h"
27
28 #include <iostream>
29
30 using namespace std;
31
32 ClassImp(AliAnalysisTaskTotEt)
33
34 //________________________________________________________________________
35 AliAnalysisTaskTotEt::AliAnalysisTaskTotEt(const char *name, Bool_t isMc) :
36         AliAnalysisTaskTransverseEnergy(name, isMc)
37         ,fRecAnalysis(0)
38         ,fMCAnalysis(0)
39 {
40     // Constructor
41
42     // select if we should use EMCal or PHOS class
43     // PHOS by default, EMCal if name string contains EMC
44     TString t(name);
45     t.ToUpper();
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()");
51         }
52
53         if (fRecoConfigFile.Length()) {
54             cout<<"Rereading AliAnalysisEtReconstructedEmcal configuration file..."<<endl;
55             gROOT->LoadMacro(fRecoConfigFile);
56             fRecAnalysis = (AliAnalysisEtReconstructedEmcal *) gInterpreter->ProcessLine("ConfigEtReconstructed()");
57         }
58     }
59     else {
60         if (fMCConfigFile.Length()) {
61             cout<<"Rereading AliAnalysisEtMonteCarloPhos configuration file..."<<endl;
62             gROOT->LoadMacro(fMCConfigFile);
63             fMCAnalysis = (AliAnalysisEtMonteCarloPhos *) gInterpreter->ProcessLine("ConfigEtMonteCarlo()");
64         }
65
66         if (fRecoConfigFile.Length()) {
67             cout<<"Rereading AliAnalysisEtReconstructedPhos configuration file..."<<endl;
68             gROOT->LoadMacro(fRecoConfigFile);
69             fRecAnalysis = (AliAnalysisEtReconstructedPhos *) gInterpreter->ProcessLine("ConfigEtReconstructed()");
70         }
71     }
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
76
77     DefineOutput(1, TList::Class());
78
79 }
80 AliAnalysisTaskTotEt::~AliAnalysisTaskTotEt() {//Destructor
81     fOutputList->Clear();
82     delete fRecAnalysis;
83     delete fMCAnalysis;
84 }
85
86 //________________________________________________________________________
87 void AliAnalysisTaskTotEt::UserCreateOutputObjects()
88 {
89     // Create histograms
90     // Called once
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);
99
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");
110     }
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");
121     }
122
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());
133         //add ITS stuff!
134     }
135     else {
136         Printf("Error: no track cuts!");
137     }
138
139 }
140
141 //________________________________________________________________________
142 void AliAnalysisTaskTotEt::UserExec(Option_t *)
143 { // execute method
144
145     fESDEvent = dynamic_cast<AliESDEvent*>(InputEvent());
146     if (!fESDEvent)
147     {
148         Printf("ERROR: Could not retrieve event");
149         return;
150     }
151
152     Int_t res = CheckPhysicsSelection(fESDEvent->GetRunNumber());
153
154     AliCentrality *cent = GetCentralityObject();
155
156     if (res == 0 && cent)
157     {
158         if (IsPhysicsSelected())
159         {
160             fRecAnalysis->SetCentralityObject(cent);
161             fRecAnalysis->AnalyseEvent(fESDEvent);
162
163             AliMCEvent* mcEvent = MCEvent();
164             if (mcEvent)
165             {
166                 fMCAnalysis->AnalyseEvent(mcEvent);
167             }
168             fHistEtRecvsEtMC->Fill(fRecAnalysis->GetTotEtAcc(), fMCAnalysis->GetTotEt());
169         }
170     }
171     // Post output data.
172     PostData(1, fOutputList);
173 }
174
175 //________________________________________________________________________
176 void AliAnalysisTaskTotEt::Terminate(Option_t *)
177 {
178     // Draw result to the screen
179     // Called once at the end of the query
180
181     fOutputList = dynamic_cast<TList*> (GetOutputData(1));
182     if (!fOutputList) {
183         printf("ERROR: Output list not available\n");
184         return;
185     }
186 }
187
188
189