]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/totEt/AliAnalysisTaskHadEt.cxx
Fixing big memory leaks from empty destructors and fixing some smaller memory leaks
[u/mrichter/AliRoot.git] / PWG4 / 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
27 #include <iostream>
28
29 using namespace std;
30
31 ClassImp(AliAnalysisTaskHadEt)
32
33
34
35 //________________________________________________________________________
36 AliAnalysisTaskHadEt::AliAnalysisTaskHadEt(const char *name, Bool_t isMc) :
37         AliAnalysisTaskTransverseEnergy(name, isMc)
38         ,fRecAnalysis(0)
39         ,fMCAnalysis(0)
40 {
41     // Constructor
42   fMCConfigFile = "ConfigHadEtMonteCarlo.C";
43   fRecoConfigFile = "ConfigHadEtReconstructed.C";
44   
45   if(fMCAnalysis) delete fMCAnalysis;
46   if(fRecAnalysis) delete fRecAnalysis;
47
48   if (fMCConfigFile.Length()) {
49     cout<<"Rereading AliAnalysisHadEtMonteCarlo configuration file..."<<endl;
50     gROOT->LoadMacro(fMCConfigFile);
51     fMCAnalysis = (AliAnalysisHadEtMonteCarlo *) gInterpreter->ProcessLine("ConfigHadEtMonteCarlo()");
52   }
53
54   if (fRecoConfigFile.Length()) {
55     cout<<"Rereading AliAnalysisHadEtReconstructed configuration file..."<<endl;
56     gROOT->LoadMacro(fRecoConfigFile);
57     fRecAnalysis = (AliAnalysisHadEtReconstructed *) gInterpreter->ProcessLine("ConfigHadEtReconstructed()");
58   }
59
60     // Define input and output slots here
61     // Input slot #0 works with a TChain
62     DefineInput(0, TChain::Class());
63     // Output slot #1 writes into a TH1 container
64
65     DefineOutput(1, TList::Class());
66 }
67 AliAnalysisTaskHadEt::~AliAnalysisTaskHadEt(){//Destructor
68   delete fRecAnalysis;
69   delete fMCAnalysis;
70 }
71
72
73 //________________________________________________________________________
74 void AliAnalysisTaskHadEt::UserCreateOutputObjects()
75 {
76     // Create histograms
77     // Called once
78   fOutputList = new TList;
79   fOutputList->SetOwner();
80   fMCAnalysis->SetHistoList(fOutputList);
81   fRecAnalysis->SetHistoList(fOutputList);
82   fMCAnalysis->CreateHistograms();
83   fRecAnalysis->CreateHistograms();
84
85   if(fRecAnalysis->DataSet() != fMCAnalysis->DataSet()){
86     cout<<"Warning: Reconstruction data set and Monte Carlo data set are not the same!  Setting data set to "<<fRecAnalysis->DataSet()<<endl;
87   }
88
89   Bool_t selectPrimaries=kTRUE;
90   if(fEsdtrackCutsITSTPC) delete fEsdtrackCutsITSTPC;
91   if(fEsdtrackCutsITS) delete fEsdtrackCutsITS;
92   if(fEsdtrackCutsTPC) delete fEsdtrackCutsTPC;
93   if(fRecAnalysis->DataSet()==2009){
94     cout<<"Setting track cuts for the 2009 p+p collisions at 900 GeV"<<endl;
95     fEsdtrackCutsITSTPC = AliESDtrackCuts::GetStandardITSTPCTrackCuts2009(selectPrimaries);
96     fEsdtrackCutsITSTPC->SetName("fEsdTrackCuts");
97     fEsdtrackCutsTPC = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
98     fEsdtrackCutsTPC->SetName("fEsdTrackCutsTPCOnly");
99     //ITS stand alone cuts - similar to 2009 cuts but with only ITS hits required
100     fEsdtrackCutsITS =  AliESDtrackCuts::GetStandardITSPureSATrackCuts2009(kTRUE,kFALSE);//we do want primaries but we do not want to require PID info
101     fEsdtrackCutsITS->SetName("fEsdTrackCutsITS");
102   }
103   if(fRecAnalysis->DataSet()==2010 || fRecAnalysis->DataSet()==20100){
104     //cout<<"Setting track cuts for the 2010 p+p collisions at 7 GeV"<<endl;
105     cout<<"Warning:  Have not set 2010 track cuts yet!!"<<endl;
106     fEsdtrackCutsITSTPC = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010(selectPrimaries);
107     fEsdtrackCutsITSTPC->SetName("fEsdTrackCuts");
108     fEsdtrackCutsTPC = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
109     fEsdtrackCutsTPC->SetName("fEsdTrackCutsTPCOnly");
110     //ITS stand alone cuts - similar to 2009 cuts but with only ITS hits required
111     fEsdtrackCutsITS =  AliESDtrackCuts::GetStandardITSPureSATrackCuts2009(kTRUE,kFALSE);//we do want primaries but we do not want to require PID info
112     fEsdtrackCutsITS->SetName("fEsdTrackCutsITS");
113   }
114
115   fOutputList->Add(fEsdtrackCutsITSTPC);
116   fOutputList->Add(fEsdtrackCutsTPC);
117   fOutputList->Add(fEsdtrackCutsITS);
118   if(fEsdtrackCutsITSTPC && fEsdtrackCutsTPC){
119     fRecAnalysis->SetITSTrackCuts( GetITSTrackCuts());
120     fMCAnalysis->SetITSTrackCuts( GetITSTrackCuts());
121     fRecAnalysis->SetTPCITSTrackCuts( GetTPCITSTrackCuts());
122     fMCAnalysis->SetTPCITSTrackCuts( GetTPCITSTrackCuts());
123     fRecAnalysis->SetTPCOnlyTrackCuts( GetTPCOnlyTrackCuts());
124     fMCAnalysis->SetTPCOnlyTrackCuts( GetTPCOnlyTrackCuts());
125     //add ITS stuff!
126   }
127   else{
128     Printf("Error: no track cuts!");
129   }
130 }
131
132 //________________________________________________________________________
133 void AliAnalysisTaskHadEt::UserExec(Option_t *)
134 { // execute method
135   fESDEvent = dynamic_cast<AliESDEvent*>(InputEvent());
136 if (!fESDEvent) {
137   Printf("ERROR: Could not retrieve event");
138   return;
139  }
140 //cout<<"AliAnalysisHadEtReconstructed 90"<<endl;
141
142 Int_t res = CheckPhysicsSelection(fESDEvent->GetRunNumber()); // Check if the physics selection is valid for this run
143
144 AliESDCentrality *cent = GetCentralityObject();
145
146 if(res == 0 && cent){
147   
148   
149   fRecAnalysis->AnalyseEvent(fESDEvent);
150
151   AliMCEvent* mcEvent = MCEvent();
152 // if (!mcEvent) {
153 //   Printf("ERROR: Could not retrieve MC event");
154 //  }
155   if (mcEvent && fESDEvent){
156       ((AliAnalysisHadEtMonteCarlo*)fMCAnalysis)->AnalyseEvent((AliVEvent*)mcEvent,(AliVEvent*)fESDEvent);
157       if(fMCAnalysis->Full()){
158         fMCAnalysis->FillSimTotEtMinusRecoTotEtFullAcceptanceTPC( fRecAnalysis->GetCorrectedTotEtFullAcceptanceTPC() );
159         fMCAnalysis->FillSimTotEtMinusRecoTotEtFullAcceptanceITS( fRecAnalysis->GetCorrectedTotEtFullAcceptanceITS() );
160         fMCAnalysis->FillSimTotEtMinusRecoTotEtFullAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedTotEtFullAcceptanceTPCNoPID() );
161         fMCAnalysis->FillSimTotEtMinusRecoTotEtFullAcceptanceITSNoPID( fRecAnalysis->GetCorrectedTotEtFullAcceptanceITSNoPID() );
162         fMCAnalysis->FillSimHadEtMinusRecoHadEtFullAcceptanceTPC( fRecAnalysis->GetCorrectedHadEtFullAcceptanceTPC() );
163         fMCAnalysis->FillSimHadEtMinusRecoHadEtFullAcceptanceITS( fRecAnalysis->GetCorrectedHadEtFullAcceptanceITS() );
164         fMCAnalysis->FillSimHadEtMinusRecoHadEtFullAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedHadEtFullAcceptanceTPCNoPID() );
165         fMCAnalysis->FillSimHadEtMinusRecoHadEtFullAcceptanceITSNoPID( fRecAnalysis->GetCorrectedHadEtFullAcceptanceITSNoPID() );
166       }
167       if(fMCAnalysis->EMCAL()){
168         fMCAnalysis->FillSimTotEtMinusRecoTotEtEMCALAcceptanceTPC( fRecAnalysis->GetCorrectedTotEtEMCALAcceptanceTPC() );
169         fMCAnalysis->FillSimTotEtMinusRecoTotEtEMCALAcceptanceITS( fRecAnalysis->GetCorrectedTotEtEMCALAcceptanceITS() );
170         fMCAnalysis->FillSimTotEtMinusRecoTotEtEMCALAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedTotEtEMCALAcceptanceTPCNoPID() );
171         fMCAnalysis->FillSimTotEtMinusRecoTotEtEMCALAcceptanceITSNoPID( fRecAnalysis->GetCorrectedTotEtEMCALAcceptanceITSNoPID() );
172         fMCAnalysis->FillSimHadEtMinusRecoHadEtEMCALAcceptanceTPC( fRecAnalysis->GetCorrectedHadEtEMCALAcceptanceTPC() );
173         fMCAnalysis->FillSimHadEtMinusRecoHadEtEMCALAcceptanceITS( fRecAnalysis->GetCorrectedHadEtEMCALAcceptanceITS() );
174         fMCAnalysis->FillSimHadEtMinusRecoHadEtEMCALAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedHadEtEMCALAcceptanceTPCNoPID() );
175         fMCAnalysis->FillSimHadEtMinusRecoHadEtEMCALAcceptanceITSNoPID( fRecAnalysis->GetCorrectedHadEtEMCALAcceptanceITSNoPID() );
176       }
177       if(fMCAnalysis->PHOS()){
178         fMCAnalysis->FillSimTotEtMinusRecoTotEtPHOSAcceptanceTPC( fRecAnalysis->GetCorrectedTotEtPHOSAcceptanceTPC() );
179         fMCAnalysis->FillSimTotEtMinusRecoTotEtPHOSAcceptanceITS( fRecAnalysis->GetCorrectedTotEtPHOSAcceptanceITS() );
180         fMCAnalysis->FillSimTotEtMinusRecoTotEtPHOSAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedTotEtPHOSAcceptanceTPCNoPID() );
181         fMCAnalysis->FillSimTotEtMinusRecoTotEtPHOSAcceptanceITSNoPID( fRecAnalysis->GetCorrectedTotEtPHOSAcceptanceITSNoPID() );
182         fMCAnalysis->FillSimHadEtMinusRecoHadEtPHOSAcceptanceTPC( fRecAnalysis->GetCorrectedHadEtPHOSAcceptanceTPC() );
183         fMCAnalysis->FillSimHadEtMinusRecoHadEtPHOSAcceptanceITS( fRecAnalysis->GetCorrectedHadEtPHOSAcceptanceITS() );
184         fMCAnalysis->FillSimHadEtMinusRecoHadEtPHOSAcceptanceTPCNoPID( fRecAnalysis->GetCorrectedHadEtPHOSAcceptanceTPCNoPID() );
185         fMCAnalysis->FillSimHadEtMinusRecoHadEtPHOSAcceptanceITSNoPID( fRecAnalysis->GetCorrectedHadEtPHOSAcceptanceITSNoPID() );
186       }
187       if(fMCAnalysis->PiKP() && fMCAnalysis->Full()){
188         fMCAnalysis->FillSimPiKPMinusRecoPiKPFullAcceptanceTPC(fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceTPC());
189         fMCAnalysis->FillSimPiKPMinusRecoPiKPFullAcceptanceITS(fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceITS());
190         fMCAnalysis->FillSimPiKPMinusRecoPiKPFullAcceptanceTPCNoPID(fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceTPCNoPID());
191         fMCAnalysis->FillSimPiKPMinusRecoPiKPFullAcceptanceITSNoPID(fRecAnalysis->GetCorrectedPiKPEtFullAcceptanceITSNoPID());
192       }
193     }
194   }
195 // Post output data.
196 PostData(1, fOutputList);
197 }
198
199 //________________________________________________________________________
200 void AliAnalysisTaskHadEt::Terminate(Option_t *)
201 {
202     // Draw result to the screen
203     // Called once at the end of the query
204
205     fOutputList = dynamic_cast<TList*> (GetOutputData(1));
206     if (!fOutputList) {
207         printf("ERROR: Output list not available\n");
208         return;
209     }
210 }
211
212