]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FEMTOSCOPY/AliFemto/AliAnalysisTaskFemto.cxx
Make the macros run for PROOF and PROOFLite
[u/mrichter/AliRoot.git] / PWG2 / FEMTOSCOPY / AliFemto / AliAnalysisTaskFemto.cxx
1 //------------------------------------------------------\r
2 // AliAnalysisTaskFemto - A task for the analysis framework\r
3 // from the FEMTOSCOPY analysis of PWG2. Creates the necessary\r
4 // connection between the ESD or AOD input and the femtoscopic\r
5 // code.\r
6 // Author: Adam Kisiel, OSU; Adam.Kisiel@cern.ch\r
7 //------------------------------------------------------\r
8 #include "TROOT.h"\r
9 #include "TChain.h"\r
10 #include "TH1.h"\r
11 #include "TCanvas.h"\r
12 #include "TSystem.h"\r
13 #include "TFile.h"\r
14 #include "TInterpreter.h"\r
15 \r
16 #include "AliAnalysisTask.h"\r
17 \r
18 #include "AliESDEvent.h"\r
19 \r
20 #include "AliFemtoAnalysis.h"\r
21 #include "AliAnalysisTaskFemto.h"\r
22 #include "AliVHeader.h"\r
23 #include "AliGenEventHeader.h"\r
24 #include "AliGenHijingEventHeader.h"\r
25 #include "AliGenCocktailEventHeader.h"\r
26 \r
27 ClassImp(AliAnalysisTaskFemto)\r
28 \r
29 // Default name for the setup macro of femto analysis  \r
30 // This function MUST be defined in the separate file !!!\r
31 // extern AliFemtoManager *ConfigFemtoAnalysis();\r
32 \r
33 //________________________________________________________________________\r
34   AliAnalysisTaskFemto::AliAnalysisTaskFemto(const char *name): \r
35     AliAnalysisTask(name,""), \r
36     fESD(0), \r
37     fAOD(0),\r
38     fStack(0),\r
39     fOutputList(0), \r
40     fReader(0x0),\r
41     fManager(0x0),\r
42     fAnalysisType(0)\r
43 {\r
44   // Constructor.\r
45   // Input slot #0 works with an Ntuple\r
46   DefineInput(0, TChain::Class());\r
47   // Output slot #0 writes into a TH1 container\r
48   DefineOutput(0, TList::Class());\r
49 }\r
50 \r
51 AliAnalysisTaskFemto::AliAnalysisTaskFemto(const AliAnalysisTaskFemto& aFemtoTask):\r
52     AliAnalysisTask(aFemtoTask), \r
53     fESD(0), \r
54     fAOD(0),\r
55     fStack(0),\r
56     fOutputList(0), \r
57     fReader(0x0),\r
58     fManager(0x0),\r
59     fAnalysisType(0)\r
60 {\r
61   // copy constructor\r
62   fESD = aFemtoTask.fESD; \r
63   fAOD = aFemtoTask.fAOD; \r
64   fStack = aFemtoTask.fStack;\r
65   fOutputList = aFemtoTask.fOutputList;   \r
66   fReader = aFemtoTask.fReader;       \r
67   fManager = aFemtoTask.fManager;      \r
68   fAnalysisType = aFemtoTask.fAnalysisType; \r
69 }\r
70 \r
71 \r
72 AliAnalysisTaskFemto& AliAnalysisTaskFemto::operator=(const AliAnalysisTaskFemto& aFemtoTask){\r
73   // assignment operator\r
74   if (this == &aFemtoTask)\r
75     return *this;\r
76 \r
77   fESD = aFemtoTask.fESD; \r
78   fAOD = aFemtoTask.fAOD; \r
79   fStack = aFemtoTask.fStack;\r
80   fOutputList = aFemtoTask.fOutputList;   \r
81   fReader = aFemtoTask.fReader;       \r
82   fManager = aFemtoTask.fManager;      \r
83   fAnalysisType = aFemtoTask.fAnalysisType; \r
84 \r
85   return *this;\r
86 }\r
87 \r
88 //________________________________________________________________________\r
89 void AliAnalysisTaskFemto::ConnectInputData(Option_t *) {\r
90   printf("   ConnectInputData %s\n", GetName());\r
91 \r
92   fESD = 0;\r
93   fAOD = 0;\r
94   fAnalysisType = 0;\r
95 \r
96   TTree* tree = dynamic_cast<TTree*> (GetInputData(0));\r
97   if (!tree) {\r
98     Printf("ERROR: Could not read chain from input slot 0");\r
99   } \r
100   else {\r
101     AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());\r
102     \r
103     if(esdH) {\r
104       cout << "Selected ESD analysis" << endl;\r
105       fAnalysisType = 1;\r
106       \r
107       if (!esdH) {\r
108         Printf("ERROR: Could not get ESDInputHandler");\r
109       } \r
110       else {\r
111         fESD = esdH->GetEvent();\r
112       }\r
113     }\r
114     else {\r
115       AliAODInputHandler *aodH = dynamic_cast<AliAODInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());\r
116       \r
117       if (!aodH) {\r
118         Printf("ERROR: Could not get AODInputHandler");\r
119       } \r
120       else {\r
121         cout << "Selected AOD analysis" << endl;\r
122         fAnalysisType = 2;\r
123 \r
124         fAOD = aodH->GetEvent();\r
125       }\r
126     }\r
127     if ((!fAOD) && (!fESD)) {\r
128       Printf("Wrong analysis type: Only ESD and AOD types are allowed!");\r
129     }\r
130   }\r
131   \r
132   \r
133 }\r
134 \r
135 //________________________________________________________________________\r
136 void AliAnalysisTaskFemto::CreateOutputObjects() {\r
137   printf("Creating Femto Analysis objects\n");\r
138 \r
139   gROOT->LoadMacro("ConfigFemtoAnalysis.C");\r
140   //  fJetFinder = (AliJetFinder*) gInterpreter->ProcessLine("ConfigJetAnalysis()");\r
141   SetFemtoManager((AliFemtoManager *) gInterpreter->ProcessLine("ConfigFemtoAnalysis()"));\r
142 \r
143   TList *tOL;\r
144   fOutputList = fManager->Analysis(0)->GetOutputList();\r
145 \r
146   for (unsigned int ian = 1; ian<fManager->AnalysisCollection()->size(); ian++) {\r
147     tOL = fManager->Analysis(ian)->GetOutputList();\r
148 \r
149     TIter nextListCf(tOL);\r
150     while (TObject *obj = nextListCf()) {\r
151       fOutputList->Add(obj);\r
152     }\r
153 \r
154     delete tOL;\r
155   }\r
156 }\r
157 \r
158 //________________________________________________________________________\r
159 void AliAnalysisTaskFemto::Exec(Option_t *) {\r
160   // Task making a femtoscopic analysis.\r
161 \r
162   if (fAnalysisType==1) {\r
163     if (!fESD) {\r
164       Printf("ERROR: fESD not available");\r
165       return;\r
166     }\r
167 \r
168     //Get MC data\r
169     AliMCEventHandler*    mctruth = (AliMCEventHandler*) \r
170       ((AliAnalysisManager::GetAnalysisManager())->GetMCtruthEventHandler());\r
171     \r
172     AliGenHijingEventHeader *hdh = 0;\r
173     if(mctruth) {\r
174       fStack = mctruth->MCEvent()->Stack();\r
175 \r
176       AliGenCocktailEventHeader *hd = dynamic_cast<AliGenCocktailEventHeader *> (mctruth->MCEvent()->GenEventHeader());\r
177       \r
178       if (hd) {\r
179         \r
180         printf ("Got MC cocktail event header %p\n", (void *) hd);\r
181         TList *lhd = hd->GetHeaders();\r
182         printf ("Got list of headers %d\n", lhd->GetEntries());\r
183         \r
184         for (int iterh=0; iterh<lhd->GetEntries(); iterh++) \r
185           {\r
186             hdh = dynamic_cast<AliGenHijingEventHeader *> (lhd->At(iterh));\r
187             printf ("HIJING header at %i is %p\n", iterh, (void *) hdh);\r
188           }\r
189       }    \r
190     }\r
191 \r
192     // Get ESD\r
193     AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());\r
194     \r
195     if (!esdH) {\r
196       Printf("ERROR: Could not get ESDInputHandler");\r
197       return;\r
198     } \r
199     else {\r
200       fESD = esdH->GetEvent();\r
201     }\r
202 \r
203     printf("Tracks in ESD: %d \n",fESD->GetNumberOfTracks());\r
204 \r
205     if (fESD->GetNumberOfTracks() >= 0) {\r
206     \r
207       if (!fReader) {\r
208         printf("ERROR: No ESD reader for ESD analysis !\n");\r
209       }\r
210       \r
211       AliFemtoEventReaderESDChain* fesdc = dynamic_cast<AliFemtoEventReaderESDChain *> (fReader);\r
212       if (fesdc)\r
213         {\r
214           // Process the event with no Kine information\r
215           fesdc->SetESDSource(fESD);\r
216           fManager->ProcessEvent();\r
217         }\r
218       AliFemtoEventReaderESDChainKine* fesdck = dynamic_cast<AliFemtoEventReaderESDChainKine *> (fReader);\r
219       if (fesdck) \r
220         {\r
221           // Process the event with Kine information\r
222           fesdck->SetESDSource(fESD);\r
223           fesdck->SetStackSource(fStack);\r
224           \r
225           fesdck->SetGenEventHeader(hdh);\r
226           fManager->ProcessEvent();\r
227         }\r
228     } \r
229 \r
230     // Post the output histogram list\r
231     PostData(0, fOutputList);\r
232   }\r
233   \r
234   if (fAnalysisType==2) {    \r
235     if (!fAOD) {\r
236       Printf("ERROR: fAOD not available");\r
237       return;\r
238     }\r
239 \r
240     // Get AOD\r
241     AliAODInputHandler *aodH = dynamic_cast<AliAODInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());\r
242       \r
243     if (!aodH) {\r
244       Printf("ERROR: Could not get AODInputHandler");\r
245       return;\r
246     } \r
247     else {\r
248 \r
249       fAOD = aodH->GetEvent();\r
250     }\r
251 \r
252     printf("Tracks in AOD: %d \n",fAOD->GetNumberOfTracks());\r
253     \r
254     if (fAOD->GetNumberOfTracks() > 0) {\r
255       if (!fReader) {\r
256         printf("ERROR: No AOD reader for AOD analysis! \n");\r
257       }\r
258       else {\r
259         AliFemtoEventReaderAODChain* faodc = dynamic_cast<AliFemtoEventReaderAODChain *> (fReader);\r
260 \r
261         if (faodc) {\r
262           // Process the event\r
263           faodc->SetAODSource(fAOD);\r
264           fManager->ProcessEvent();\r
265         }\r
266       }\r
267     } \r
268 \r
269     // Post the output histogram list\r
270     PostData(0, fOutputList);\r
271   }\r
272 }      \r
273 \r
274 //________________________________________________________________________\r
275 void AliAnalysisTaskFemto::Terminate(Option_t *) {\r
276   // Do the final processing\r
277   if (fManager) {\r
278     fManager->Finish();\r
279   }\r
280 }\r
281 //________________________________________________________________________\r
282 void AliAnalysisTaskFemto:: FinishTaskOutput() {\r
283   // Do the final processing\r
284   if (fManager) {\r
285     fManager->Finish();\r
286   }\r
287 }\r
288 //________________________________________________________________________\r
289 void AliAnalysisTaskFemto::SetFemtoReaderESD(AliFemtoEventReaderESDChain *aReader)\r
290 {\r
291   printf("Selectring Femto reader for ESD\n");\r
292   fReader = aReader;\r
293 }\r
294 //________________________________________________________________________\r
295 void AliAnalysisTaskFemto::SetFemtoReaderESDKine(AliFemtoEventReaderESDChainKine *aReader)\r
296 {\r
297   printf("Selectring Femto reader for ESD with Kinematics information\n");\r
298   fReader = aReader;\r
299 }\r
300 //________________________________________________________________________\r
301 void AliAnalysisTaskFemto::SetFemtoReaderAOD(AliFemtoEventReaderAODChain *aReader)\r
302 {\r
303   printf("Selecting Femto reader for AOD\n");\r
304   fReader = aReader;\r
305 }\r
306 //________________________________________________________________________\r
307 void AliAnalysisTaskFemto::SetFemtoManager(AliFemtoManager *aManager)\r
308 {\r
309   fManager = aManager;\r
310   printf("Got reader %p\n", (void *) aManager->EventReader());\r
311   AliFemtoEventReaderESDChain     *tReaderESDChain     = dynamic_cast<AliFemtoEventReaderESDChain *> (aManager->EventReader());\r
312   AliFemtoEventReaderESDChainKine *tReaderESDChainKine = dynamic_cast<AliFemtoEventReaderESDChainKine *> (aManager->EventReader());\r
313   AliFemtoEventReaderAODChain     *tReaderAODChain     = dynamic_cast<AliFemtoEventReaderAODChain *> (aManager->EventReader());\r
314 \r
315   if ((!tReaderESDChain) && (!tReaderESDChainKine) && (!tReaderAODChain)) {\r
316     printf("No AliFemto event reader created. Will not run femto analysis.\n");\r
317     return;\r
318   }\r
319   if (tReaderESDChain) SetFemtoReaderESD(tReaderESDChain);\r
320   if (tReaderESDChainKine) SetFemtoReaderESDKine(tReaderESDChainKine);\r
321   if (tReaderAODChain) SetFemtoReaderAOD(tReaderAODChain);\r
322 }\r
323 \r