]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FLOW/AliFlowTasks/AliAnalysisTaskFlowEvent.cxx
use a seperate config file to set the cuts and options for the flow analysis
[u/mrichter/AliRoot.git] / PWG2 / FLOW / AliFlowTasks / AliAnalysisTaskFlowEvent.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" //needed as include
20 #include "TList.h"
21
22 // ALICE Analysis Framework
23 class AliAnalysisTask;
24 #include "AliAnalysisManager.h"
25
26 // ESD interface
27 #include "AliESDEvent.h"
28 #include "AliESDInputHandler.h"
29
30 // AOD interface
31 #include "AliAODEvent.h"
32 #include "AliAODInputHandler.h"
33
34 // Monte Carlo Event
35 #include "AliMCEventHandler.h"
36 #include "AliMCEvent.h"
37
38 // ALICE Correction Framework
39 #include "AliCFManager.h"
40
41 // Interface to Event generators to get Reaction Plane Angle
42 #include "AliGenCocktailEventHeader.h"
43 #include "AliGenHijingEventHeader.h"
44 #include "AliGenGeVSimEventHeader.h"
45
46 // Interface to make the Flow Event Simple used in the flow analysis methods
47 #include "AliFlowEventSimpleMaker.h"
48
49 #include "AliAnalysisTaskFlowEvent.h"
50
51 // AliAnalysisTaskFlowEvent:
52 //
53 // analysis task for filling the flow event
54 // from MCEvent, ESD, AOD ....
55 // and put it in an output stream so it can 
56 // be used by the various flow analysis methods 
57 // for cuts the correction framework is used
58 // which also outputs QA histrograms to view
59 // the effects of the cuts
60
61
62 ClassImp(AliAnalysisTaskFlowEvent)
63   
64 extern void ConfigFlowAnalysis(AliAnalysisTaskFlowEvent*); 
65
66 //________________________________________________________________________
67 AliAnalysisTaskFlowEvent::AliAnalysisTaskFlowEvent(const char *name, Bool_t on) : 
68   AliAnalysisTask(name, ""), 
69 //  fOutputFile(NULL),
70   fESD(NULL),
71   fAOD(NULL),
72   fEventMaker(NULL),
73   fAnalysisType("ESD"),
74   fCFManager1(NULL),
75   fCFManager2(NULL),
76   fQAInt(NULL),
77   fQADiff(NULL),
78   fQA(on)
79 {
80   // Constructor
81   cout<<"AliAnalysisTaskFlowEvent::AliAnalysisTaskFlowEvent(const char *name)"<<endl;
82
83   // Define input and output slots here
84   // Input slot #0 works with a TChain
85   DefineInput(0, TChain::Class());
86   // Define here the flow event output
87   DefineOutput(0, AliFlowEventSimple::Class());  
88   if(on) {
89     DefineOutput(1, TList::Class());
90     DefineOutput(2, TList::Class()); }  
91   // and for testing open an output file
92   //  fOutputFile = new TFile("FlowEvents.root","RECREATE");
93
94 }
95
96 //________________________________________________________________________
97 AliAnalysisTaskFlowEvent::AliAnalysisTaskFlowEvent() : 
98   //  fOutputFile(NULL),
99   fESD(NULL),
100   fAOD(NULL),
101   fEventMaker(NULL),
102   fAnalysisType("ESD"),
103   fCFManager1(NULL),
104   fCFManager2(NULL),
105   fQAInt(NULL),
106   fQADiff(NULL),
107   fQA(kFALSE)
108 {
109   // Constructor
110   cout<<"AliAnalysisTaskFlowEvent::AliAnalysisTaskFlowEvent()"<<endl;
111 }
112
113 //________________________________________________________________________
114 AliAnalysisTaskFlowEvent::~AliAnalysisTaskFlowEvent()
115 {
116   //
117   // Destructor
118   //
119
120   // objects in the output list are deleted 
121   // by the TSelector dtor (I hope)
122
123 }
124
125 //________________________________________________________________________
126 void AliAnalysisTaskFlowEvent::ConnectInputData(Option_t *) 
127 {
128   // Connect ESD or AOD here
129   // Called once
130   cout<<"AliAnalysisTaskFlowEvent::ConnectInputData(Option_t *)"<<endl;
131
132   TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
133   if (!tree) {
134     Printf("ERROR: Could not read chain from input slot 0");
135   } 
136   else {
137     if (fAnalysisType == "ESD" || fAnalysisType == "ESDMC0" || fAnalysisType == "ESDMC1" || fAnalysisType == "MC") {
138       AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
139       if (!esdH) {
140         Printf("ERROR: Could not get ESDInputHandler");
141       } 
142       else {fESD = esdH->GetEvent();}
143     }
144     else if (fAnalysisType == "AOD") {
145       AliAODInputHandler *aodH = dynamic_cast<AliAODInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
146       if (!aodH) {
147         Printf("ERROR: Could not get AODInputHandler");
148       }
149       else { fAOD = aodH->GetEvent();}
150     }
151     else {
152       Printf("!!!!!Wrong analysis type: Only ESD, ESDMC0, ESDMC1, AOD and MC types are allowed!");
153       exit(1);
154     }
155   }
156 }
157
158 //________________________________________________________________________
159 void AliAnalysisTaskFlowEvent::CreateOutputObjects() 
160 {
161   // Called at every worker node to initialize
162   cout<<"AliAnalysisTaskFlowEvent::CreateOutputObjects()"<<endl;
163
164   if (!(fAnalysisType == "AOD" || fAnalysisType == "ESD" || fAnalysisType == "ESDMC0"  || fAnalysisType == "ESDMC1" || fAnalysisType == "MC")) {
165     cout<<"WRONG ANALYSIS TYPE! only ESD, ESDMC0, ESDMC1, AOD and MC are allowed."<<endl;
166     exit(1);
167   }
168
169   ConfigFlowAnalysis(this); 
170
171   // Flow Event maker
172   fEventMaker = new AliFlowEventSimpleMaker();
173 }
174
175 //________________________________________________________________________
176 void AliAnalysisTaskFlowEvent::Exec(Option_t *) 
177 {
178   // Main loop
179   // Called for each event
180   AliFlowEventSimple* fEvent = NULL;
181   Double_t fRP = 0.; // the monte carlo reaction plane angle
182   AliMCEvent* mcEvent = NULL;
183   // See if we can get Monte Carlo Information and if so get the reaction plane
184
185   AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
186   if (eventHandler) {
187     mcEvent = eventHandler->MCEvent();
188     if (mcEvent) {
189       if (!strcmp(mcEvent-> GenEventHeader()->GetName(),"Cocktail Header")) { //returns 0 if matches
190         AliGenCocktailEventHeader *headerC = dynamic_cast<AliGenCocktailEventHeader *> (mcEvent-> GenEventHeader()); 
191         if (headerC) {
192           TList *lhd = headerC->GetHeaders();
193           if (lhd) {
194             AliGenHijingEventHeader *hdh = dynamic_cast<AliGenHijingEventHeader *> (lhd->At(0)); 
195             if (hdh) {
196               fRP = hdh->ReactionPlaneAngle();
197               //cout<<"The reactionPlane from Hijing is: "<< fRP <<endl;
198             }
199           }
200         }
201         //else { cout<<"headerC is NULL"<<endl; }
202       }
203       else if (!strcmp(mcEvent-> GenEventHeader()->GetName(),"GeVSim header")) { //returns 0 if matches
204         AliGenGeVSimEventHeader* headerG = (AliGenGeVSimEventHeader*)(mcEvent->GenEventHeader());
205         if (headerG) {
206           fRP = headerG->GetEventPlane();
207           //cout<<"The reactionPlane from GeVSim is: "<< fRP <<endl;
208         }
209         //else { cout<<"headerG is NULL"<<endl; }
210       }
211     }
212     //else {cout<<"No MC event!"<<endl; }
213   }
214   //else {cout<<"No eventHandler!"<<endl; }
215
216   // set the value of the monte carlo event plane for the flow event
217   fEventMaker->SetMCReactionPlaneAngle(fRP);
218
219   // Fill the FlowEventSimple for MC input          
220   if (fAnalysisType == "MC") {
221     if (!fCFManager1) {cout << "ERROR: No pointer to correction framework cuts! " << endl; return; }
222     if (!fCFManager2) {cout << "ERROR: No pointer to correction framework cuts! " << endl; return; }
223
224     // Process MC truth, therefore we receive the AliAnalysisManager and ask it for the AliMCEventHandler
225     // This handler can return the current MC event
226     if (!mcEvent) { Printf("ERROR: Could not retrieve MC event"); return;}
227
228     fCFManager1->SetEventInfo(mcEvent);
229     fCFManager2->SetEventInfo(mcEvent);
230
231     // analysis 
232     Printf("Number of MC particles: %d", mcEvent->GetNumberOfTracks());
233     fEvent = fEventMaker->FillTracks(mcEvent,fCFManager1,fCFManager2);
234     // here we have the fEvent and want to make it available as an output stream
235     // so no delete fEvent;
236   }
237   // Fill the FlowEventSimple for ESD input  
238   else if (fAnalysisType == "ESD") {
239     if (!fCFManager1) {cout << "ERROR: No pointer to correction framework cuts! " << endl; return; }
240     if (!fCFManager2) {cout << "ERROR: No pointer to correction framework cuts! " << endl; return; }
241
242     if (!fESD) { Printf("ERROR: fESD not available"); return;}
243     Printf("There are %d tracks in this event", fESD->GetNumberOfTracks());
244     
245     // analysis 
246     fEvent = fEventMaker->FillTracks(fESD,fCFManager1,fCFManager2);
247   }
248   // Fill the FlowEventSimple for ESD input combined with MC info  
249   else if (fAnalysisType == "ESDMC0" || fAnalysisType == "ESDMC1" ) {
250     if (!fCFManager1) {cout << "ERROR: No pointer to correction framework cuts! " << endl; return; }
251     if (!fCFManager2) {cout << "ERROR: No pointer to correction framework cuts! " << endl; return; }
252     if (!fESD) { Printf("ERROR: fESD not available"); return;}
253     Printf("There are %d tracks in this event", fESD->GetNumberOfTracks());
254     
255     if (!mcEvent) {Printf("ERROR: Could not retrieve MC event"); return;}
256
257     fCFManager1->SetEventInfo(mcEvent);
258     fCFManager2->SetEventInfo(mcEvent);
259
260
261     if (fAnalysisType == "ESDMC0") { 
262       fEvent = fEventMaker->FillTracks(fESD, mcEvent, fCFManager1, fCFManager2, 0); //0 = kine from ESD, 1 = kine from MC
263     } else if (fAnalysisType == "ESDMC1") {
264       fEvent = fEventMaker->FillTracks(fESD, mcEvent, fCFManager1, fCFManager2, 1); //0 = kine from ESD, 1 = kine from MC
265     }
266   }
267   // Fill the FlowEventSimple for AOD input  
268   else if (fAnalysisType == "AOD") {
269     if (!fAOD) {Printf("ERROR: fAOD not available"); return;}
270     Printf("There are %d tracks in this event", fAOD->GetNumberOfTracks());
271
272     // analysis 
273     //For the moment don't use CF //AliFlowEventSimple* fEvent = fEventMaker->FillTracks(fAOD,fCFManager1,fCFManager2);
274     fEvent = fEventMaker->FillTracks(fAOD);
275   }
276
277   //fListHistos->Print();
278   //  fOutputFile->WriteObject(fEvent,"myFlowEventSimple");     
279   PostData(0,fEvent);
280   if (fQA) {
281     PostData(1,fQAInt);
282     PostData(2,fQADiff); }
283
284
285 //________________________________________________________________________
286 void AliAnalysisTaskFlowEvent::Terminate(Option_t *) 
287 {
288   // Called once at the end of the query -- do not call in case of CAF
289
290 }