]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG2/FLOW/AliAnalysisTaskLYZEventPlane.cxx
last one
[u/mrichter/AliRoot.git] / PWG2 / FLOW / AliAnalysisTaskLYZEventPlane.cxx
CommitLineData
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"
20
21class AliAnalysisTask;
22#include "AliAnalysisManager.h"
23
24#include "AliESDEvent.h"
25#include "AliESDInputHandler.h"
26
27#include "AliAODEvent.h"
28#include "AliAODInputHandler.h"
29
30#include "AliMCEventHandler.h"
31#include "AliMCEvent.h"
32
33#include "AliAnalysisTaskLYZEventPlane.h"
34#include "AliFlowEventSimpleMaker.h"
35#include "AliFlowLYZEventPlane.h"
36#include "AliFlowAnalysisWithLYZEventPlane.h"
37
38// AliAnalysisTaskLYZEventPlane:
39//
40// analysis task for Lee Yang Zeros Event Plane
41//
42// Author: Naomi van der Kolk (kolk@nikhef.nl)
43
44ClassImp(AliAnalysisTaskLYZEventPlane)
45
46//________________________________________________________________________
47AliAnalysisTaskLYZEventPlane::AliAnalysisTaskLYZEventPlane(const char *name) :
48 AliAnalysisTask(name, ""),
d97a8fc5 49 fESD(NULL),
50 fAOD(NULL),
f1d945a1 51 fAnalysisType("ESD"),
d97a8fc5 52 fLyzEp(NULL),
53 fLyz(NULL),
54 fEventMaker(NULL),
55 fFirstRunFile(NULL),
56 fSecondRunFile(NULL)
f1d945a1 57{
58 // Constructor
59 cout<<"AliAnalysisTaskLYZEventPlane::AliAnalysisTaskLYZEventPlane(const char *name)"<<endl;
60
61 // Define input and output slots here
62 // Input slot #0 works with a TChain
63 DefineInput(0, TChain::Class());
64 DefineInput(1, TList::Class());
65 DefineInput(2, TList::Class());
66 // Output slot #0 writes into a TList container
67 DefineOutput(0, TList::Class());
68}
69
70//________________________________________________________________________
71void AliAnalysisTaskLYZEventPlane::ConnectInputData(Option_t *)
72{
73 // Connect ESD or AOD here
74 // Called once
75 cout<<"AliAnalysisTaskLYZEventPlane::ConnectInputData(Option_t *)"<<endl;
76
77 TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
78 if (!tree) {
79 Printf("ERROR: Could not read chain from input slot 0");
80 } else {
81 // Disable all branches and enable only the needed ones
82 if (fAnalysisType == "MC") {
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 }
0b7f49e9 94 else if (fAnalysisType == "ESD" || fAnalysisType == "ESDMC0" || fAnalysisType == "ESDMC1" ) {
f1d945a1 95 tree->SetBranchStatus("*", kFALSE);
96 tree->SetBranchStatus("Tracks.*", kTRUE);
97
98 AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
99
100 if (!esdH) {
101 Printf("ERROR: Could not get ESDInputHandler");
102 } else
103 fESD = esdH->GetEvent();
104 }
105 else if (fAnalysisType == "AOD") {
106 AliAODInputHandler *aodH = dynamic_cast<AliAODInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
107
108 if (!aodH) {
109 Printf("ERROR: Could not get AODInputHandler");
110 }
111 else {
112 fAOD = aodH->GetEvent();
113 }
114 }
115 else {
0b7f49e9 116 Printf("Wrong analysis type: Only ESD, ESDMC0, ESDMC1, AOD and MC types are allowed!");
f1d945a1 117 exit(1);
118
119 }
120 }
121}
122
123//________________________________________________________________________
124void AliAnalysisTaskLYZEventPlane::CreateOutputObjects()
125{
126 // Called once
127 cout<<"AliAnalysisTaskLYZEventPlane::CreateOutputObjects()"<<endl;
128
0b7f49e9 129 if (!(fAnalysisType == "AOD" || fAnalysisType == "ESD" || fAnalysisType == "ESDMC0" || fAnalysisType == "ESDMC1" || fAnalysisType == "MC")) {
130 cout<<"WRONG ANALYSIS TYPE! only ESD, , ESDMC0, ESDMC1, AOD and MC are allowed."<<endl;
f1d945a1 131 exit(1);
132 }
133
134 //event maker
135 fEventMaker = new AliFlowEventSimpleMaker();
136 //lee yang zeros event plane
137 fLyzEp = new AliFlowLYZEventPlane() ;
138 //Analyser
139 fLyz = new AliFlowAnalysisWithLYZEventPlane() ;
140
141 //output file
fa92c1a0 142 TString outputName = "outputFromLYZEventPlaneAnalysis" ;
143 outputName += fAnalysisType.Data() ;
144 outputName += ".root" ;
145 fLyz->SetOutFileName( outputName.Data() );
f1d945a1 146
147 // Get data from input slot 1 and 2
148 fFirstRunFile = (TFile*)GetInputData(1);
149 cerr<<"fFirstRunFile ("<<fFirstRunFile<<")"<<endl;
150 if (fFirstRunFile) cerr<<"fFirstRunFile -> IsOpen() = "<<fFirstRunFile -> IsOpen()<<endl;
151 fSecondRunFile = (TFile*)GetInputData(2);
152 cerr<<"fSecondRunFile ("<<fSecondRunFile<<")"<<endl;
153 if (fSecondRunFile) cerr<<"fSecondRunFile -> IsOpen() = "<<fSecondRunFile -> IsOpen()<<endl;
154
155 fLyzEp -> SetFirstRunFile(fFirstRunFile);
156 fLyzEp -> SetSecondRunFile(fSecondRunFile);
157 fLyz -> SetFirstRunFile(fFirstRunFile);
158 fLyz -> SetSecondRunFile(fSecondRunFile);
159
160 fLyzEp-> Init();
161 fLyz-> Init();
162
163
164}
165
166//________________________________________________________________________
167void AliAnalysisTaskLYZEventPlane::Exec(Option_t *)
168{
169 // Main loop
170 // Called for each event
171 if (fAnalysisType == "MC") {
172 // Process MC truth, therefore we receive the AliAnalysisManager and ask it for the AliMCEventHandler
173 // This handler can return the current MC event
174
175 AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
176 if (!eventHandler) {
177 Printf("ERROR: Could not retrieve MC event handler");
178 return;
179 }
180
181 AliMCEvent* mcEvent = eventHandler->MCEvent();
182 if (!mcEvent) {
183 Printf("ERROR: Could not retrieve MC event");
184 return;
185 }
186
187 Printf("MC particles: %d", mcEvent->GetNumberOfTracks());
188
189 //lee yang zeros analysis
190 AliFlowEventSimple* fEvent = fEventMaker->FillTracks(mcEvent);
191 fLyz->Make(fEvent,fLyzEp);
192
193 delete fEvent;
194 }
195 else if (fAnalysisType == "ESD") {
196 if (!fESD) {
197 Printf("ERROR: fESD not available");
198 return;
199 }
200 Printf("There are %d tracks in this event", fESD->GetNumberOfTracks());
201
202 //lee yang zeros analysis
203 AliFlowEventSimple* fEvent = fEventMaker->FillTracks(fESD);
204 fLyz->Make(fEvent,fLyzEp);
205 delete fEvent;
206 }
0b7f49e9 207 else if (fAnalysisType == "ESDMC0") {
208 if (!fESD) {
209 Printf("ERROR: fESD not available");
210 return;
211 }
212 Printf("There are %d tracks in this event", fESD->GetNumberOfTracks());
213
214 AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
215 if (!eventHandler) {
216 Printf("ERROR: Could not retrieve MC event handler");
217 return;
218 }
219
220 AliMCEvent* mcEvent = eventHandler->MCEvent();
221 if (!mcEvent) {
222 Printf("ERROR: Could not retrieve MC event");
223 return;
224 }
225
226 //lee yang zeros analysis
227 AliFlowEventSimple* fEvent = fEventMaker->FillTracks(fESD,mcEvent,0); //0 = kine from ESD, 1 = kine from MC
228 fLyz->Make(fEvent,fLyzEp);
229 delete fEvent;
230 //delete mcEvent;
231 }
232 else if (fAnalysisType == "ESDMC1") {
233 if (!fESD) {
234 Printf("ERROR: fESD not available");
235 return;
236 }
237 Printf("There are %d tracks in this event", fESD->GetNumberOfTracks());
238
239 AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
240 if (!eventHandler) {
241 Printf("ERROR: Could not retrieve MC event handler");
242 return;
243 }
244
245 AliMCEvent* mcEvent = eventHandler->MCEvent();
246 if (!mcEvent) {
247 Printf("ERROR: Could not retrieve MC event");
248 return;
249 }
250
251 //lee yang zeros analysis
252 AliFlowEventSimple* fEvent = fEventMaker->FillTracks(fESD,mcEvent,1); //0 = kine from ESD, 1 = kine from MC
253 fLyz->Make(fEvent,fLyzEp);
254 delete fEvent;
255 //delete mcEvent;
256 }
f1d945a1 257 else if (fAnalysisType == "AOD") {
258 if (!fAOD) {
259 Printf("ERROR: fAOD not available");
260 return;
261 }
262 Printf("There are %d tracks in this event", fAOD->GetNumberOfTracks());
263
264 //lee yang zeros analysis
265 AliFlowEventSimple* fEvent = fEventMaker->FillTracks(fAOD);
266 fLyz->Make(fEvent,fLyzEp);
267 delete fEvent;
268 }
269
270 PostData(0,fLyz->GetOutFile());
271}
272
273//________________________________________________________________________
274void AliAnalysisTaskLYZEventPlane::Terminate(Option_t *)
275{
276 // Called once at the end of the query
277 fLyz->Finish();
278
279 delete fLyz;
280 delete fLyzEp;
281 delete fEventMaker;
282}
af795c87 283
284