]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FLOW/AliFlowTasks/AliAnalysisTaskLYZEventPlane.cxx
bug in StepManager fixed, thanks to Andreas
[u/mrichter/AliRoot.git] / PWG2 / FLOW / AliFlowTasks / AliAnalysisTaskLYZEventPlane.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 "TProfile.h"
20 #include "TFile.h"
21 #include "TList.h"
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 "AliCFManager.h"
36
37 #include "AliAnalysisTaskLYZEventPlane.h"
38 #include "AliFlowEventSimpleMaker.h"
39 #include "AliFlowCommonHist.h"
40 #include "AliFlowCommonHistResults.h"
41 #include "AliFlowLYZEventPlane.h"
42 #include "AliFlowAnalysisWithLYZEventPlane.h"
43
44 // AliAnalysisTaskLYZEventPlane:
45 //
46 // analysis task for Lee Yang Zeros Event Plane
47 //
48 // Author: Naomi van der Kolk (kolk@nikhef.nl)
49
50 ClassImp(AliAnalysisTaskLYZEventPlane)
51
52 //________________________________________________________________________
53 AliAnalysisTaskLYZEventPlane::AliAnalysisTaskLYZEventPlane(const char *name, Bool_t on) : 
54   AliAnalysisTask(name, ""), 
55   fESD(NULL), 
56   fAOD(NULL),
57   fAnalysisType("ESD"), 
58   fLyzEp(NULL),
59   fLyz(NULL),
60   fEventMaker(NULL),
61   fCFManager1(NULL),
62   fCFManager2(NULL),
63   fListHistos(NULL),
64   fSecondRunFile(NULL),
65   fQAInt(NULL),
66   fQADiff(NULL),
67   fQA(on)
68 {
69   // Constructor
70   cout<<"AliAnalysisTaskLYZEventPlane::AliAnalysisTaskLYZEventPlane(const char *name)"<<endl;
71
72   // Define input and output slots here
73   // Input slot #0 works with a TChain
74   DefineInput(0, TChain::Class());
75   DefineInput(1, TList::Class());
76   // Output slot #0 writes into a TList container
77   DefineOutput(0, TList::Class());
78   if(on) {
79     DefineOutput(1, TList::Class());
80     DefineOutput(2, TList::Class()); } 
81 }
82
83 //________________________________________________________________________
84 AliAnalysisTaskLYZEventPlane::AliAnalysisTaskLYZEventPlane() : 
85   fESD(NULL), 
86   fAOD(NULL),
87   fAnalysisType("ESD"), 
88   fLyzEp(NULL),
89   fLyz(NULL),
90   fEventMaker(NULL),
91   fCFManager1(NULL),
92   fCFManager2(NULL),
93   fListHistos(NULL),
94   fSecondRunFile(NULL),
95   fQAInt(NULL),
96   fQADiff(NULL),
97   fQA(kFALSE)
98 {
99   // Constructor
100   cout<<"AliAnalysisTaskLYZEventPlane::AliAnalysisTaskLYZEventPlane()"<<endl;
101 }
102
103
104 //________________________________________________________________________
105 AliAnalysisTaskLYZEventPlane::~AliAnalysisTaskLYZEventPlane() 
106 {
107   //destructor
108
109 }
110
111 //________________________________________________________________________
112 void AliAnalysisTaskLYZEventPlane::ConnectInputData(Option_t *) 
113 {
114   // Connect ESD or AOD here
115   // Called once
116   cout<<"AliAnalysisTaskLYZEventPlane::ConnectInputData(Option_t *)"<<endl;
117
118   TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
119   if (!tree) {
120     Printf("ERROR: Could not read chain from input slot 0");
121   } else {
122     // Disable all branches and enable only the needed ones
123     if (fAnalysisType == "MC") {
124       // we want to process only MC
125       tree->SetBranchStatus("*", kFALSE);
126
127       AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
128
129       if (!esdH) {
130         Printf("ERROR: Could not get ESDInputHandler");
131       } else {
132         fESD = esdH->GetEvent();
133       }
134     }
135     else if (fAnalysisType == "ESD" || fAnalysisType == "ESDMC0" || fAnalysisType == "ESDMC1" ) {
136       tree->SetBranchStatus("*", kFALSE);
137       tree->SetBranchStatus("Tracks.*", kTRUE);
138
139       AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
140
141       if (!esdH) {
142         Printf("ERROR: Could not get ESDInputHandler");
143       } else 
144         fESD = esdH->GetEvent();
145     }
146     else if (fAnalysisType == "AOD") {
147       AliAODInputHandler *aodH = dynamic_cast<AliAODInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
148
149       if (!aodH) {
150         Printf("ERROR: Could not get AODInputHandler");
151       }
152       else {
153         fAOD = aodH->GetEvent();
154       }
155     }
156     else {
157       Printf("Wrong analysis type: Only ESD, ESDMC0, ESDMC1, AOD and MC types are allowed!");
158       exit(1);
159
160     }
161   }  
162 }
163
164 //________________________________________________________________________
165 void AliAnalysisTaskLYZEventPlane::CreateOutputObjects() 
166 {
167   // Called once
168   cout<<"AliAnalysisTaskLYZEventPlane::CreateOutputObjects()"<<endl;
169
170   if (!(fAnalysisType == "AOD" || fAnalysisType == "ESD"  || fAnalysisType == "ESDMC0"  || fAnalysisType == "ESDMC1" || fAnalysisType == "MC")) {
171     cout<<"WRONG ANALYSIS TYPE! only ESD, , ESDMC0, ESDMC1, AOD and MC are allowed."<<endl;
172     exit(1);
173   }
174
175   //event maker
176   fEventMaker = new AliFlowEventSimpleMaker();
177   //lee yang zeros event plane
178   fLyzEp = new AliFlowLYZEventPlane() ;
179   //Analyser
180   fLyz = new AliFlowAnalysisWithLYZEventPlane() ;
181      
182   // Get data from input slot
183   TList* pSecondRunList = (TList*)GetInputData(1);
184   if (pSecondRunList) {
185     fLyzEp -> SetSecondRunList(pSecondRunList);
186     fLyz -> SetSecondRunList(pSecondRunList);
187   } else { cout<<"No Second run List!"<<endl; exit(0); }
188
189   fLyzEp-> Init();
190   fLyz-> Init();
191
192   if (fLyz->GetHistList()) {
193     fListHistos = fLyz->GetHistList();
194     fListHistos->Print();
195   }
196   else { cout<<"ERROR: Could not retrieve histogram list"<<endl;}
197
198
199 }
200
201 //________________________________________________________________________
202 void AliAnalysisTaskLYZEventPlane::Exec(Option_t *) 
203 {
204   // Main loop
205   // Called for each event
206
207   if (fAnalysisType == "MC") {
208     // Process MC truth, therefore we receive the AliAnalysisManager and ask it for the AliMCEventHandler
209     // This handler can return the current MC event
210
211     AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
212     if (!eventHandler) {
213       Printf("ERROR: Could not retrieve MC event handler");
214       return;
215     }
216
217     AliMCEvent* mcEvent = eventHandler->MCEvent();
218     if (!mcEvent) {
219       Printf("ERROR: Could not retrieve MC event");
220       return;
221     }
222
223     fCFManager1->SetEventInfo(mcEvent);
224     fCFManager2->SetEventInfo(mcEvent);
225
226     Printf("MC particles: %d", mcEvent->GetNumberOfTracks());
227
228     //lee yang zeros analysis 
229     AliFlowEventSimple* fEvent = fEventMaker->FillTracks(mcEvent,fCFManager1,fCFManager2);
230     fLyz->Make(fEvent,fLyzEp);
231
232     delete fEvent;
233   }
234   else if (fAnalysisType == "ESD") {
235     if (!fESD) {
236       Printf("ERROR: fESD not available");
237       return;
238     }
239     Printf("There are %d tracks in this event", fESD->GetNumberOfTracks());
240     
241     //lee yang zeros analysis 
242     AliFlowEventSimple* fEvent = fEventMaker->FillTracks(fESD,fCFManager1,fCFManager2);
243     fLyz->Make(fEvent,fLyzEp);
244     delete fEvent;
245   }
246   else if (fAnalysisType == "ESDMC0" || fAnalysisType == "ESDMC1") {
247     if (!fESD) {
248       Printf("ERROR: fESD not available");
249       return;
250     }
251     Printf("There are %d tracks in this event", fESD->GetNumberOfTracks());
252     
253     AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
254     if (!eventHandler) {
255       Printf("ERROR: Could not retrieve MC event handler");
256       return;
257     }
258
259     AliMCEvent* mcEvent = eventHandler->MCEvent();
260     if (!mcEvent) {
261       Printf("ERROR: Could not retrieve MC event");
262       return;
263     }
264
265     fCFManager1->SetEventInfo(mcEvent);
266     fCFManager2->SetEventInfo(mcEvent);
267
268     //lee yang zeros analysis 
269     AliFlowEventSimple* fEvent = NULL;
270     if (fAnalysisType == "ESDMC0") {
271       fEvent = fEventMaker->FillTracks(fESD, mcEvent,fCFManager1, fCFManager2,0);  //0 = kine from ESD, 1 = kine from MC
272     } else if (fAnalysisType == "ESDMC1") {
273       fEvent = fEventMaker->FillTracks(fESD, mcEvent,fCFManager1, fCFManager2,0);  //0 = kine from ESD, 1 = kine from MC
274     }
275
276     fLyz->Make(fEvent,fLyzEp);
277     delete fEvent;
278     //delete mcEvent;
279   }
280   
281   else if (fAnalysisType == "AOD") {
282     if (!fAOD) {
283       Printf("ERROR: fAOD not available");
284       return;
285     }
286     Printf("There are %d tracks in this event", fAOD->GetNumberOfTracks());
287
288     //lee yang zeros analysis 
289     //for the moment don't use CF
290     //AliFlowEventSimple* fEvent = fEventMaker->FillTracks(fAOD,fCFManager1,fCFManager2);
291     AliFlowEventSimple* fEvent = fEventMaker->FillTracks(fAOD);
292     fLyz->Make(fEvent,fLyzEp);
293     delete fEvent;
294   }
295   
296   PostData(0,fListHistos);
297   if (fQA) {
298     PostData(1,fQAInt);
299     PostData(2,fQADiff); }
300 }      
301
302 //________________________________________________________________________
303 void AliAnalysisTaskLYZEventPlane::Terminate(Option_t *) 
304 {
305   // Called once at the end of the query
306   AliFlowAnalysisWithLYZEventPlane* fLyzTerm = new AliFlowAnalysisWithLYZEventPlane() ;
307   fListHistos = (TList*)GetOutputData(0);
308   //cout << "histogram list in Terminate" << endl;
309   if (fListHistos) {
310     //Get the common histograms from the output list
311     AliFlowCommonHist *pCommonHist = dynamic_cast<AliFlowCommonHist*> 
312       (fListHistos->FindObject("AliFlowCommonHistLYZEP"));
313     AliFlowCommonHistResults *pCommonHistResults = dynamic_cast<AliFlowCommonHistResults*> 
314       (fListHistos->FindObject("AliFlowCommonHistResultsLYZEP"));
315
316     TProfile* pHistProR0theta = dynamic_cast<TProfile*> 
317       (fListHistos->FindObject("First_FlowPro_r0theta_LYZ"));
318
319     TProfile* pHistProFlow = dynamic_cast<TProfile*> 
320       (fListHistos->FindObject("FlowPro_VPt_LYZEP"));
321
322     TH1F* pHistQsumforChi = dynamic_cast<TH1F*> 
323       (fListHistos->FindObject("Flow_QsumforChi_LYZEP"));
324
325     if (pCommonHist && pCommonHistResults && pHistProR0theta &&
326         pHistProFlow && pHistQsumforChi ) {
327     fLyzTerm->SetCommonHists(pCommonHist);
328     fLyzTerm->SetCommonHistsRes(pCommonHistResults);
329     fLyzTerm->SetFirstr0theta(pHistProR0theta);
330     fLyzTerm->SetHistProFlow(pHistProFlow);
331     fLyzTerm->SetHistQsumforChi(pHistQsumforChi);
332     fLyzTerm->Finish();
333     PostData(0,fListHistos);
334     } else { 
335       cout<<"WARNING: Histograms needed to run Finish() are not accessable!"<<endl; 
336     }
337
338     fListHistos->Print(); 
339   } else { cout << "histogram list pointer is empty" << endl;}
340
341   cout<<".....finished LYZ EventPlane"<<endl;  
342 }
343
344