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