1 /*************************************************************************
2 * Copyright(c) 1998-2008, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 #include "Riostream.h" //needed as include
21 class AliAnalysisTask;
22 #include "AliAnalysisManager.h"
24 #include "AliESDEvent.h"
25 #include "AliESDInputHandler.h"
27 #include "AliAODEvent.h"
28 #include "AliAODInputHandler.h"
30 #include "AliMCEventHandler.h"
31 #include "AliMCEvent.h"
33 #include "AliAnalysisTaskLYZEventPlane.h"
34 #include "AliFlowEventSimpleMaker.h"
35 #include "AliFlowLYZEventPlane.h"
36 #include "AliFlowAnalysisWithLYZEventPlane.h"
38 // AliAnalysisTaskLYZEventPlane:
40 // analysis task for Lee Yang Zeros Event Plane
42 // Author: Naomi van der Kolk (kolk@nikhef.nl)
44 ClassImp(AliAnalysisTaskLYZEventPlane)
46 //________________________________________________________________________
47 AliAnalysisTaskLYZEventPlane::AliAnalysisTaskLYZEventPlane(const char *name) :
48 AliAnalysisTask(name, ""),
59 cout<<"AliAnalysisTaskLYZEventPlane::AliAnalysisTaskLYZEventPlane(const char *name)"<<endl;
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());
70 //________________________________________________________________________
71 void AliAnalysisTaskLYZEventPlane::ConnectInputData(Option_t *)
73 // Connect ESD or AOD here
75 cout<<"AliAnalysisTaskLYZEventPlane::ConnectInputData(Option_t *)"<<endl;
77 TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
79 Printf("ERROR: Could not read chain from input slot 0");
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);
86 AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
89 Printf("ERROR: Could not get ESDInputHandler");
91 fESD = esdH->GetEvent();
94 else if (fAnalysisType == "ESD" || fAnalysisType == "ESDMC0" || fAnalysisType == "ESDMC1" ) {
95 tree->SetBranchStatus("*", kFALSE);
96 tree->SetBranchStatus("Tracks.*", kTRUE);
98 AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
101 Printf("ERROR: Could not get ESDInputHandler");
103 fESD = esdH->GetEvent();
105 else if (fAnalysisType == "AOD") {
106 AliAODInputHandler *aodH = dynamic_cast<AliAODInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
109 Printf("ERROR: Could not get AODInputHandler");
112 fAOD = aodH->GetEvent();
116 Printf("Wrong analysis type: Only ESD, ESDMC0, ESDMC1, AOD and MC types are allowed!");
123 //________________________________________________________________________
124 void AliAnalysisTaskLYZEventPlane::CreateOutputObjects()
127 cout<<"AliAnalysisTaskLYZEventPlane::CreateOutputObjects()"<<endl;
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;
135 fEventMaker = new AliFlowEventSimpleMaker();
136 //lee yang zeros event plane
137 fLyzEp = new AliFlowLYZEventPlane() ;
139 fLyz = new AliFlowAnalysisWithLYZEventPlane() ;
142 TString outputName = "outputFromLYZEventPlaneAnalysis" ;
143 outputName += fAnalysisType.Data() ;
144 outputName += ".root" ;
145 fLyz->SetOutFileName( outputName.Data() );
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;
155 fLyzEp -> SetFirstRunFile(fFirstRunFile);
156 fLyzEp -> SetSecondRunFile(fSecondRunFile);
157 fLyz -> SetFirstRunFile(fFirstRunFile);
158 fLyz -> SetSecondRunFile(fSecondRunFile);
166 //________________________________________________________________________
167 void AliAnalysisTaskLYZEventPlane::Exec(Option_t *)
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
175 AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
177 Printf("ERROR: Could not retrieve MC event handler");
181 AliMCEvent* mcEvent = eventHandler->MCEvent();
183 Printf("ERROR: Could not retrieve MC event");
187 Printf("MC particles: %d", mcEvent->GetNumberOfTracks());
189 //lee yang zeros analysis
190 AliFlowEventSimple* fEvent = fEventMaker->FillTracks(mcEvent);
191 fLyz->Make(fEvent,fLyzEp);
195 else if (fAnalysisType == "ESD") {
197 Printf("ERROR: fESD not available");
200 Printf("There are %d tracks in this event", fESD->GetNumberOfTracks());
202 //lee yang zeros analysis
203 AliFlowEventSimple* fEvent = fEventMaker->FillTracks(fESD);
204 fLyz->Make(fEvent,fLyzEp);
207 else if (fAnalysisType == "ESDMC0") {
209 Printf("ERROR: fESD not available");
212 Printf("There are %d tracks in this event", fESD->GetNumberOfTracks());
214 AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
216 Printf("ERROR: Could not retrieve MC event handler");
220 AliMCEvent* mcEvent = eventHandler->MCEvent();
222 Printf("ERROR: Could not retrieve MC event");
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);
232 else if (fAnalysisType == "ESDMC1") {
234 Printf("ERROR: fESD not available");
237 Printf("There are %d tracks in this event", fESD->GetNumberOfTracks());
239 AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
241 Printf("ERROR: Could not retrieve MC event handler");
245 AliMCEvent* mcEvent = eventHandler->MCEvent();
247 Printf("ERROR: Could not retrieve MC event");
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);
257 else if (fAnalysisType == "AOD") {
259 Printf("ERROR: fAOD not available");
262 Printf("There are %d tracks in this event", fAOD->GetNumberOfTracks());
264 //lee yang zeros analysis
265 AliFlowEventSimple* fEvent = fEventMaker->FillTracks(fAOD);
266 fLyz->Make(fEvent,fLyzEp);
270 PostData(0,fLyz->GetOutFile());
273 //________________________________________________________________________
274 void AliAnalysisTaskLYZEventPlane::Terminate(Option_t *)
276 // Called once at the end of the query