]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FLOW/AliAnalysisTaskLeeYangZeros.cxx
62d1207b2e94968a9b508108937affcb27fac4dc
[u/mrichter/AliRoot.git] / PWG2 / FLOW / AliAnalysisTaskLeeYangZeros.cxx
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
21
22 class AliAnalysisTask;
23 #include "AliAnalysisManager.h"
24
25 #include "AliESDEvent.h"
26 #include "AliESDInputHandler.h"
27
28 #include "AliAODEvent.h"
29 #include "AliAODInputHandler.h"
30
31 #include "AliMCEventHandler.h"
32 #include "AliMCEvent.h"
33
34 #include "AliAnalysisTaskLeeYangZeros.h"
35 #include "AliFlowEventSimpleMaker.h"
36 #include "AliFlowAnalysisWithLeeYangZeros.h"
37
38 // AliAnalysisTaskLeeYangZeros:
39 // analysis task for Lee Yang Zeros method
40 // Author: Naomi van der Kolk (kolk@nikhef.nl)
41
42 ClassImp(AliAnalysisTaskLeeYangZeros)
43
44 //________________________________________________________________________
45 AliAnalysisTaskLeeYangZeros::AliAnalysisTaskLeeYangZeros(const char *name, Bool_t firstrun) : 
46   AliAnalysisTask(name, ""), 
47   fESD(0),
48   fAOD(0),
49   fAnalysisType("ESD"), 
50   fLyz(0),
51   fEventMaker(0),
52   fFirstRunFile(0),
53   fFirstRunLYZ(firstrun), //set boolean for firstrun to initial value
54   fUseSumLYZ(kTRUE)       //set boolean for use sum to initial value
55 {
56   // Constructor
57   cout<<"AliAnalysisTaskLeeYangZeros::AliAnalysisTaskLeeYangZeros(const char *name)"<<endl;
58
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());  
65 }
66
67 //________________________________________________________________________
68 void AliAnalysisTaskLeeYangZeros::ConnectInputData(Option_t *) 
69 {
70   // Connect ESD or AOD here
71   // Called once
72   cout<<"AliAnalysisTaskLeeYangZeros::ConnectInputData(Option_t *)"<<endl;
73
74   TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
75   if (!tree) {
76     Printf("ERROR: Could not read chain from input slot 0");
77   } 
78   else {
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);
83
84       AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
85
86       if (!esdH) {
87         Printf("ERROR: Could not get ESDInputHandler");
88       } else {
89         fESD = esdH->GetEvent();
90       }
91     }
92     else if (fAnalysisType == "ESD" || fAnalysisType == "ESDMC0" || fAnalysisType == "ESDMC1"  ) {
93       tree->SetBranchStatus("*", kFALSE);
94       tree->SetBranchStatus("Tracks.*", kTRUE);
95
96       AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
97
98       if (!esdH) {
99         Printf("ERROR: Could not get ESDInputHandler");
100       } else
101         fESD = esdH->GetEvent();
102     }
103     else if (fAnalysisType == "AOD") {
104       AliAODInputHandler *aodH = dynamic_cast<AliAODInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
105
106       if (!aodH) {
107         Printf("ERROR: Could not get AODInputHandler");
108       }
109       else {
110         fAOD = aodH->GetEvent();
111       }
112     }
113     else {
114       Printf("Wrong analysis type: Only ESD, ESDMC0, ESDMC1, AOD and MC types are allowed!");
115
116     }
117   }
118 }
119
120 //________________________________________________________________________
121 void AliAnalysisTaskLeeYangZeros::CreateOutputObjects() 
122 {
123   // Called once
124   cout<<"AliAnalysisTaskLeeYangZeros::CreateOutputObjects()"<<endl;
125
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;
128     exit(1);
129   }
130
131   //event maker
132   fEventMaker = new AliFlowEventSimpleMaker();
133   //Analyser
134   fLyz = new AliFlowAnalysisWithLeeYangZeros() ;
135    
136   fLyz -> SetFirstRun(GetFirstRunLYZ());   //set first run true or false
137   fLyz -> SetUseSum(GetUseSumLYZ());       //set use sum true or false
138
139   //output file
140   TString outputName = "outputFromLeeYangZerosAnalysis" ;
141   outputName += fAnalysisType.Data() ;
142   if (fFirstRunLYZ) {
143     outputName += "_firstrun.root" ;
144   } else {
145     outputName += "_secondrun.root" ;
146   }
147   fLyz->SetHistFileName( outputName.Data() );
148   
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;
154     
155     fLyz -> SetFirstRunFile(fFirstRunFile);
156   }
157   
158   fLyz-> Init();
159
160 }
161
162 //________________________________________________________________________
163 void AliAnalysisTaskLeeYangZeros::Exec(Option_t *) 
164 {
165   // Main loop
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
170
171     AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
172     if (!eventHandler) {
173       Printf("ERROR: Could not retrieve MC event handler");
174       return;
175     }
176
177     AliMCEvent* mcEvent = eventHandler->MCEvent();
178     if (!mcEvent) {
179       Printf("ERROR: Could not retrieve MC event");
180       return;
181     }
182
183     Printf("MC particles: %d", mcEvent->GetNumberOfTracks());
184
185     //lee yang zeros analysis 
186     AliFlowEventSimple* fEvent = fEventMaker->FillTracks(mcEvent);
187     fLyz->Make(fEvent);
188     delete fEvent;
189   }
190   else if (fAnalysisType == "ESD") {
191     if (!fESD) {
192       Printf("ERROR: fESD not available");
193       return;
194     }
195     Printf("There are %d tracks in this event", fESD->GetNumberOfTracks());
196     
197     //lee yang zeros analysis 
198     AliFlowEventSimple* fEvent = fEventMaker->FillTracks(fESD);
199     fLyz->Make(fEvent);
200     delete fEvent;
201   }
202   else if (fAnalysisType == "ESDMC0") {
203     if (!fESD) {
204       Printf("ERROR: fESD not available");
205       return;
206     }
207     Printf("There are %d tracks in this event", fESD->GetNumberOfTracks());
208     
209     AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
210     if (!eventHandler) {
211       Printf("ERROR: Could not retrieve MC event handler");
212       return;
213     }
214
215     AliMCEvent* mcEvent = eventHandler->MCEvent();
216     if (!mcEvent) {
217       Printf("ERROR: Could not retrieve MC event");
218       return;
219     }
220
221     //lee yang zeros analysis 
222     AliFlowEventSimple* fEvent = fEventMaker->FillTracks(fESD,mcEvent,0); //0 = kine from ESD, 1 = kine from MC
223     fLyz->Make(fEvent);
224     delete fEvent;
225     //delete mcEvent;
226   }
227   else if (fAnalysisType == "ESDMC1") {
228     if (!fESD) {
229       Printf("ERROR: fESD not available");
230       return;
231     }
232     Printf("There are %d tracks in this event", fESD->GetNumberOfTracks());
233     
234     AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
235     if (!eventHandler) {
236       Printf("ERROR: Could not retrieve MC event handler");
237       return;
238     }
239
240     AliMCEvent* mcEvent = eventHandler->MCEvent();
241     if (!mcEvent) {
242       Printf("ERROR: Could not retrieve MC event");
243       return;
244     }
245
246     //lee yang zeros analysis 
247     AliFlowEventSimple* fEvent = fEventMaker->FillTracks(fESD,mcEvent,1); //0 = kine from ESD, 1 = kine from MC
248     fLyz->Make(fEvent);
249     delete fEvent;
250     //delete mcEvent;
251   }
252   else if (fAnalysisType == "AOD") {
253     if (!fAOD) {
254       Printf("ERROR: fAOD not available");
255       return;
256     }
257     Printf("There are %d tracks in this event", fAOD->GetNumberOfTracks());
258
259     //lee yang zeros analysis 
260     AliFlowEventSimple* fEvent = fEventMaker->FillTracks(fAOD);
261     fLyz->Make(fEvent);
262     delete fEvent;
263   }
264   
265 }      
266
267 //________________________________________________________________________
268 void AliAnalysisTaskLeeYangZeros::Terminate(Option_t *) 
269 {
270   // Called once at the end of the query
271   if (GetNinputs() == 2) { 
272     cerr<<"fFirstRunFile -> IsOpen() = "<<fFirstRunFile -> IsOpen()<<endl;
273   }
274   cerr<<"fLyz->GetHistFile() -> IsOpen() = "<<fLyz->GetHistFile() -> IsOpen()<<endl;
275
276   fLyz->Finish();
277
278   PostData(0,fLyz->GetHistFile());
279
280   delete fLyz;
281   delete fEventMaker;
282
283   cout<<".....finished"<<endl;
284 }