]>
Commit | Line | Data |
---|---|---|
f1d945a1 | 1 | /************************************************************************* |
2 | * Copyright(c) 1998-2008, ALICE Experiment at CERN, All rights reserved. * | |
3 | * * | |
4 | * Author: The ALICE Off-line Project. * | |
5 | * Contributors are mentioned in the code where appropriate. * | |
6 | * * | |
7 | * Permission to use, copy, modify and distribute this software and its * | |
8 | * documentation strictly for non-commercial purposes is hereby granted * | |
9 | * without fee, provided that the above copyright notice appears in all * | |
10 | * copies and that both the copyright notice and this permission notice * | |
11 | * appear in the supporting documentation. The authors make no claims * | |
12 | * about the suitability of this software for any purpose. It is * | |
13 | * provided "as is" without express or implied warranty. * | |
14 | **************************************************************************/ | |
15 | ||
16 | #include "Riostream.h" //needed as include | |
17 | #include "TChain.h" | |
18 | #include "TTree.h" | |
19 | #include "TFile.h" //needed as include | |
20 | #include "TList.h" | |
21 | ||
22 | ||
23 | class AliAnalysisTask; | |
24 | #include "AliAnalysisManager.h" | |
25 | ||
26 | #include "AliESDEvent.h" | |
27 | #include "AliESDInputHandler.h" | |
28 | ||
29 | #include "AliAODEvent.h" | |
30 | #include "AliAODInputHandler.h" | |
31 | ||
32 | #include "AliMCEventHandler.h" | |
33 | #include "AliMCEvent.h" | |
34 | ||
35 | #include "AliGenCocktailEventHeader.h" | |
36 | #include "AliGenHijingEventHeader.h" | |
37 | ||
38 | #include "AliAnalysisTaskMCEventPlane.h" | |
39 | #include "AliFlowEventSimpleMaker.h" | |
40 | #include "AliFlowAnalysisWithMCEventPlane.h" | |
41 | ||
42 | // AliAnalysisTaskMCEventPlane: | |
43 | // | |
44 | // analysis task for Monte Carlo Event Plane | |
45 | // | |
46 | // Author: Naomi van der Kolk (kolk@nikhef.nl) | |
47 | ||
48 | ClassImp(AliAnalysisTaskMCEventPlane) | |
49 | ||
50 | //________________________________________________________________________ | |
51 | AliAnalysisTaskMCEventPlane::AliAnalysisTaskMCEventPlane(const char *name) : | |
52 | AliAnalysisTask(name, ""), | |
53 | fESD(0), | |
54 | fAOD(0), | |
55 | fMc(0), | |
56 | fEventMaker(0), | |
57 | fAnalysisType("MC") | |
58 | { | |
59 | // Constructor | |
60 | cout<<"AliAnalysisTaskMCEventPlane::AliAnalysisTaskMCEventPlane(const char *name)"<<endl; | |
61 | ||
62 | // Define input and output slots here | |
63 | // Input slot #0 works with a TChain | |
64 | DefineInput(0, TChain::Class()); | |
65 | // Output slot #0 writes into a TList container | |
66 | DefineOutput(0, TList::Class()); | |
67 | } | |
68 | ||
69 | //________________________________________________________________________ | |
70 | void AliAnalysisTaskMCEventPlane::ConnectInputData(Option_t *) | |
71 | { | |
72 | // Connect ESD or AOD here | |
73 | // Called once | |
74 | cout<<"AliAnalysisTaskMCEventPlane::ConnectInputData(Option_t *)"<<endl; | |
75 | ||
76 | TTree* tree = dynamic_cast<TTree*> (GetInputData(0)); | |
77 | if (!tree) { | |
78 | Printf("ERROR: Could not read chain from input slot 0"); | |
79 | } else { | |
80 | // Disable all branches and enable only the needed ones | |
81 | if (fAnalysisType == "MC") { | |
82 | cout<<"!!!!!reading MC kinematics only"<<endl; | |
83 | // we want to process only MC | |
84 | tree->SetBranchStatus("*", kFALSE); | |
85 | ||
86 | AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()); | |
87 | ||
88 | if (!esdH) { | |
89 | Printf("ERROR: Could not get ESDInputHandler"); | |
90 | } else { | |
91 | fESD = esdH->GetEvent(); | |
92 | } | |
93 | } | |
94 | else if (fAnalysisType == "ESD") { | |
95 | cout<<"!!!!!reading the ESD only"<<endl; | |
96 | tree->SetBranchStatus("*", kFALSE); | |
97 | tree->SetBranchStatus("Tracks.*", kTRUE); | |
98 | ||
99 | AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()); | |
100 | ||
101 | if (!esdH) { | |
102 | Printf("ERROR: Could not get ESDInputHandler"); | |
103 | } else | |
104 | fESD = esdH->GetEvent(); | |
105 | } | |
106 | else if (fAnalysisType == "AOD") { | |
107 | cout<<"!!!!!reading the AOD only"<<endl; | |
108 | AliAODInputHandler *aodH = dynamic_cast<AliAODInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()); | |
109 | ||
110 | if (!aodH) { | |
111 | Printf("ERROR: Could not get AODInputHandler"); | |
112 | } | |
113 | else { | |
114 | fAOD = aodH->GetEvent(); | |
115 | } | |
116 | } | |
117 | else { | |
118 | Printf("!!!!!Wrong analysis type: Only ESD, AOD and MC types are allowed!"); | |
119 | exit(1); | |
120 | ||
121 | } | |
122 | } | |
123 | } | |
124 | ||
125 | //________________________________________________________________________ | |
126 | void AliAnalysisTaskMCEventPlane::CreateOutputObjects() | |
127 | { | |
128 | // Called once | |
129 | cout<<"AliAnalysisTaskMCEventPlane::CreateOutputObjects()"<<endl; | |
130 | ||
131 | if (!(fAnalysisType == "AOD" || fAnalysisType == "ESD" || fAnalysisType == "MC")) { | |
132 | cout<<"WRONG ANALYSIS TYPE! only ESD, AOD and MC are allowed."<<endl; | |
133 | exit(1); | |
134 | } | |
135 | ||
136 | //event maker | |
137 | fEventMaker = new AliFlowEventSimpleMaker(); | |
138 | //Analyser | |
139 | fMc = new AliFlowAnalysisWithMCEventPlane() ; | |
140 | ||
141 | //output file | |
142 | TString fName = "outputFromMCEventPlaneAnalysis" ; | |
143 | fName += fAnalysisType.Data() ; | |
144 | fName += ".root" ; | |
145 | fMc->SetHistFileName( fName.Data() ); | |
146 | ||
147 | fMc-> Init(); | |
148 | ||
149 | ||
150 | } | |
151 | ||
152 | //________________________________________________________________________ | |
153 | void AliAnalysisTaskMCEventPlane::Exec(Option_t *) | |
154 | { | |
155 | // Main loop | |
156 | // Called for each event | |
157 | ||
158 | ||
159 | // Process MC truth, therefore we receive the AliAnalysisManager and ask it for the AliMCEventHandler | |
160 | // This handler can return the current MC event | |
161 | ||
162 | AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler()); | |
163 | if (!eventHandler) { | |
164 | Printf("ERROR: Could not retrieve MC event handler"); | |
165 | return; | |
166 | } | |
167 | ||
168 | AliMCEvent* mcEvent = eventHandler->MCEvent(); | |
169 | if (!mcEvent) { | |
170 | Printf("ERROR: Could not retrieve MC event"); | |
171 | return; | |
172 | } | |
173 | ||
174 | Printf("MC particles: %d", mcEvent->GetNumberOfTracks()); | |
175 | ||
176 | AliGenCocktailEventHeader *header = dynamic_cast<AliGenCocktailEventHeader *> (mcEvent-> GenEventHeader()); | |
177 | if (!header) { | |
178 | Printf("ERROR: Could not retrieve AliGenCocktailEventHeader"); | |
179 | return; | |
180 | } | |
181 | ||
182 | TList *lhd = header->GetHeaders(); | |
183 | if (!lhd) { | |
184 | Printf("ERROR: Could not retrieve List of headers"); | |
185 | return; | |
186 | } | |
187 | ||
188 | AliGenHijingEventHeader *hdh = dynamic_cast<AliGenHijingEventHeader *> (lhd->At(0)); | |
189 | if (!hdh) { | |
190 | Printf("ERROR: Could not retrieve AliGenHijingEventHeader"); | |
191 | return; | |
192 | } | |
193 | ||
194 | Double_t fRP = hdh->ReactionPlaneAngle(); | |
195 | cout<<"The reactionPlane is "<<hdh->ReactionPlaneAngle()<<endl; | |
196 | ||
197 | if (fAnalysisType == "MC") { | |
198 | // analysis | |
199 | AliFlowEventSimple* fEvent = fEventMaker->FillTracks(mcEvent); | |
200 | fMc->Make(fEvent,fRP); | |
201 | ||
202 | delete fEvent; | |
203 | } | |
204 | else if (fAnalysisType == "ESD") { | |
205 | if (!fESD) { | |
206 | Printf("ERROR: fESD not available"); | |
207 | return; | |
208 | } | |
209 | Printf("There are %d tracks in this event", fESD->GetNumberOfTracks()); | |
210 | ||
211 | // analysis | |
212 | AliFlowEventSimple* fEvent = fEventMaker->FillTracks(fESD); | |
213 | fMc->Make(fEvent,fRP); | |
214 | delete fEvent; | |
215 | } | |
216 | else if (fAnalysisType == "AOD") { | |
217 | if (!fAOD) { | |
218 | Printf("ERROR: fAOD not available"); | |
219 | return; | |
220 | } | |
221 | Printf("There are %d tracks in this event", fAOD->GetNumberOfTracks()); | |
222 | ||
223 | // analysis | |
224 | AliFlowEventSimple* fEvent = fEventMaker->FillTracks(fAOD); | |
225 | fMc->Make(fEvent,fRP); | |
226 | delete fEvent; | |
227 | } | |
228 | ||
229 | } | |
230 | ||
231 | //________________________________________________________________________ | |
232 | void AliAnalysisTaskMCEventPlane::Terminate(Option_t *) | |
233 | { | |
234 | // Called once at the end of the query | |
235 | fMc->Finish(); | |
236 | PostData(0,fMc->GetHistFile()); | |
237 | ||
238 | delete fMc; | |
239 | delete fEventMaker; | |
240 | } |