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
22 class AliAnalysisTask;
23 #include "AliAnalysisManager.h"
25 #include "AliESDEvent.h"
26 #include "AliESDInputHandler.h"
28 #include "AliAODEvent.h"
29 #include "AliAODInputHandler.h"
31 #include "AliMCEventHandler.h"
32 #include "AliMCEvent.h"
34 #include "AliAnalysisTaskLeeYangZeros.h"
35 #include "AliFlowEventSimpleMaker.h"
36 #include "AliFlowAnalysisWithLeeYangZeros.h"
38 // AliAnalysisTaskLeeYangZeros:
39 // analysis task for Lee Yang Zeros method
40 // Author: Naomi van der Kolk (kolk@nikhef.nl)
42 ClassImp(AliAnalysisTaskLeeYangZeros)
44 //________________________________________________________________________
45 AliAnalysisTaskLeeYangZeros::AliAnalysisTaskLeeYangZeros(const char *name, Bool_t firstrun) :
46 AliAnalysisTask(name, ""),
53 fFirstRunLYZ(firstrun), //set boolean for firstrun to initial value
54 fUseSumLYZ(kTRUE) //set boolean for use sum to initial value
57 cout<<"AliAnalysisTaskLeeYangZeros::AliAnalysisTaskLeeYangZeros(const char *name)"<<endl;
59 // Define input and output slots here
60 // Input slot #0 works with a TChain
61 DefineInput(0, TChain::Class());
62 if (!firstrun) DefineInput(1, TList::Class()); //for second loop
63 // Output slot #0 writes into a TList container
64 DefineOutput(0, TList::Class());
67 //________________________________________________________________________
68 void AliAnalysisTaskLeeYangZeros::ConnectInputData(Option_t *)
70 // Connect ESD or AOD here
72 cout<<"AliAnalysisTaskLeeYangZeros::ConnectInputData(Option_t *)"<<endl;
74 TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
76 Printf("ERROR: Could not read chain from input slot 0");
79 // Disable all branches and enable only the needed ones
80 if (fAnalysisType == "MC") {
81 // we want to process only MC
82 tree->SetBranchStatus("*", kFALSE);
84 AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
87 Printf("ERROR: Could not get ESDInputHandler");
89 fESD = esdH->GetEvent();
92 else if (fAnalysisType == "ESD" || fAnalysisType == "ESDMC0" || fAnalysisType == "ESDMC1" ) {
93 tree->SetBranchStatus("*", kFALSE);
94 tree->SetBranchStatus("Tracks.*", kTRUE);
96 AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
99 Printf("ERROR: Could not get ESDInputHandler");
101 fESD = esdH->GetEvent();
103 else if (fAnalysisType == "AOD") {
104 AliAODInputHandler *aodH = dynamic_cast<AliAODInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
107 Printf("ERROR: Could not get AODInputHandler");
110 fAOD = aodH->GetEvent();
114 Printf("Wrong analysis type: Only ESD, ESDMC0, ESDMC1, AOD and MC types are allowed!");
120 //________________________________________________________________________
121 void AliAnalysisTaskLeeYangZeros::CreateOutputObjects()
124 cout<<"AliAnalysisTaskLeeYangZeros::CreateOutputObjects()"<<endl;
126 if (!(fAnalysisType == "AOD" || fAnalysisType == "ESD" || fAnalysisType == "ESDMC0" || fAnalysisType == "ESDMC1" || fAnalysisType == "MC")) {
127 cout<<"WRONG ANALYSIS TYPE! only ESD, ESDMC0, ESDMC1, AOD and MC are allowed."<<endl;
132 fEventMaker = new AliFlowEventSimpleMaker();
134 fLyz = new AliFlowAnalysisWithLeeYangZeros() ;
136 fLyz -> SetFirstRun(GetFirstRunLYZ()); //set first run true or false
137 fLyz -> SetUseSum(GetUseSumLYZ()); //set use sum true or false
140 TString outputName = "outputFromLeeYangZerosAnalysis" ;
141 outputName += fAnalysisType.Data() ;
143 outputName += "_firstrun.root" ;
145 outputName += "_secondrun.root" ;
147 fLyz->SetHistFileName( outputName.Data() );
149 // Get data from input slot 1
150 if (GetNinputs() == 2) { //if there are two input slots
151 fFirstRunFile = (TFile*)GetInputData(1);
152 cerr<<"fFirstRunFile ("<<fFirstRunFile<<")"<<endl;
153 if (fFirstRunFile) cerr<<"fFirstRunFile -> IsOpen() = "<<fFirstRunFile -> IsOpen()<<endl;
155 fLyz -> SetFirstRunFile(fFirstRunFile);
162 //________________________________________________________________________
163 void AliAnalysisTaskLeeYangZeros::Exec(Option_t *)
166 // Called for each event
167 if (fAnalysisType == "MC") {
168 // Process MC truth, therefore we receive the AliAnalysisManager and ask it for the AliMCEventHandler
169 // This handler can return the current MC event
171 AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
173 Printf("ERROR: Could not retrieve MC event handler");
177 AliMCEvent* mcEvent = eventHandler->MCEvent();
179 Printf("ERROR: Could not retrieve MC event");
183 Printf("MC particles: %d", mcEvent->GetNumberOfTracks());
185 //lee yang zeros analysis
186 AliFlowEventSimple* fEvent = fEventMaker->FillTracks(mcEvent);
190 else if (fAnalysisType == "ESD") {
192 Printf("ERROR: fESD not available");
195 Printf("There are %d tracks in this event", fESD->GetNumberOfTracks());
197 //lee yang zeros analysis
198 AliFlowEventSimple* fEvent = fEventMaker->FillTracks(fESD);
202 else if (fAnalysisType == "ESDMC0") {
204 Printf("ERROR: fESD not available");
207 Printf("There are %d tracks in this event", fESD->GetNumberOfTracks());
209 AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
211 Printf("ERROR: Could not retrieve MC event handler");
215 AliMCEvent* mcEvent = eventHandler->MCEvent();
217 Printf("ERROR: Could not retrieve MC event");
221 //lee yang zeros analysis
222 AliFlowEventSimple* fEvent = fEventMaker->FillTracks(fESD,mcEvent,0); //0 = kine from ESD, 1 = kine from MC
227 else if (fAnalysisType == "ESDMC1") {
229 Printf("ERROR: fESD not available");
232 Printf("There are %d tracks in this event", fESD->GetNumberOfTracks());
234 AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
236 Printf("ERROR: Could not retrieve MC event handler");
240 AliMCEvent* mcEvent = eventHandler->MCEvent();
242 Printf("ERROR: Could not retrieve MC event");
246 //lee yang zeros analysis
247 AliFlowEventSimple* fEvent = fEventMaker->FillTracks(fESD,mcEvent,1); //0 = kine from ESD, 1 = kine from MC
252 else if (fAnalysisType == "AOD") {
254 Printf("ERROR: fAOD not available");
257 Printf("There are %d tracks in this event", fAOD->GetNumberOfTracks());
259 //lee yang zeros analysis
260 AliFlowEventSimple* fEvent = fEventMaker->FillTracks(fAOD);
267 //________________________________________________________________________
268 void AliAnalysisTaskLeeYangZeros::Terminate(Option_t *)
270 // Called once at the end of the query
271 if (GetNinputs() == 2) {
272 cerr<<"fFirstRunFile -> IsOpen() = "<<fFirstRunFile -> IsOpen()<<endl;
274 cerr<<"fLyz->GetHistFile() -> IsOpen() = "<<fLyz->GetHistFile() -> IsOpen()<<endl;
278 PostData(0,fLyz->GetHistFile());
283 cout<<".....finished"<<endl;