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