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