]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FEMTOSCOPY/AliFemto/AliAnalysisTaskFemto.cxx
Add possibility to run in ESD-to-AOD train. Fix Coverity 10087
[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, const char *aConfigMacro="ConfigFemtoAnalysis.C"): \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     fConfigMacro(0)\r
44 {\r
45   // Constructor.\r
46   // Input slot #0 works with an Ntuple\r
47   DefineInput(0, TChain::Class());\r
48   // Output slot #0 writes into a TH1 container\r
49   DefineOutput(0, TList::Class());\r
50   fConfigMacro = (char *) malloc(sizeof(char) * strlen(aConfigMacro));\r
51   strcpy(fConfigMacro, aConfigMacro);\r
52 }\r
53 \r
54 AliAnalysisTaskFemto::AliAnalysisTaskFemto(const AliAnalysisTaskFemto& aFemtoTask):\r
55     AliAnalysisTask(aFemtoTask), \r
56     fESD(0), \r
57     fAOD(0),\r
58     fStack(0),\r
59     fOutputList(0), \r
60     fReader(0x0),\r
61     fManager(0x0),\r
62     fAnalysisType(0),\r
63     fConfigMacro(0)\r
64 {\r
65   // copy constructor\r
66   fESD = aFemtoTask.fESD; \r
67   fAOD = aFemtoTask.fAOD; \r
68   fStack = aFemtoTask.fStack;\r
69   fOutputList = aFemtoTask.fOutputList;   \r
70   fReader = aFemtoTask.fReader;       \r
71   fManager = aFemtoTask.fManager;      \r
72   fAnalysisType = aFemtoTask.fAnalysisType; \r
73   fConfigMacro = (char *) malloc(sizeof(char) * strlen(aFemtoTask.fConfigMacro));\r
74   strcpy(fConfigMacro, aFemtoTask.fConfigMacro);\r
75 }\r
76 \r
77 \r
78 AliAnalysisTaskFemto& AliAnalysisTaskFemto::operator=(const AliAnalysisTaskFemto& aFemtoTask){\r
79   // assignment operator\r
80   if (this == &aFemtoTask)\r
81     return *this;\r
82 \r
83   fESD = aFemtoTask.fESD; \r
84   fAOD = aFemtoTask.fAOD; \r
85   fStack = aFemtoTask.fStack;\r
86   fOutputList = aFemtoTask.fOutputList;   \r
87   fReader = aFemtoTask.fReader;       \r
88   fManager = aFemtoTask.fManager;      \r
89   fAnalysisType = aFemtoTask.fAnalysisType; \r
90   if (fConfigMacro) free(fConfigMacro);\r
91   fConfigMacro = (char *) malloc(sizeof(char) * strlen(aFemtoTask.fConfigMacro));\r
92   strcpy(fConfigMacro, aFemtoTask.fConfigMacro);\r
93 \r
94   return *this;\r
95 }\r
96 \r
97 AliAnalysisTaskFemto::~AliAnalysisTaskFemto() \r
98 {\r
99   if (fConfigMacro) free(fConfigMacro);\r
100 }\r
101 \r
102 \r
103 //________________________________________________________________________\r
104 void AliAnalysisTaskFemto::ConnectInputData(Option_t *) {\r
105   AliInfo(Form("   ConnectInputData %s\n", GetName()));\r
106 \r
107   fESD = 0;\r
108   fAOD = 0;\r
109   fAnalysisType = 0;\r
110 \r
111   TTree* tree = dynamic_cast<TTree*> (GetInputData(0));\r
112   if (!tree) {\r
113     AliWarning("Could not read chain from input slot 0");\r
114     return;\r
115   } \r
116 \r
117   if ((dynamic_cast<AliFemtoEventReaderESDChain *> (fReader)) ||\r
118       (dynamic_cast<AliFemtoEventReaderESDChainKine *> (fReader))) {\r
119     AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());\r
120     \r
121     if(esdH) {\r
122       AliInfo("Selected ESD analysis");\r
123       fAnalysisType = 1;\r
124       \r
125 //       if (!esdH) {\r
126 //      AliWarning("Could not get ESDInputHandler");\r
127 //       } \r
128 //       else {\r
129         fESD = esdH->GetEvent();\r
130 //       }\r
131     }\r
132   }\r
133   \r
134   if (dynamic_cast<AliFemtoEventReaderAODChain *> (fReader)) {\r
135     AliAODInputHandler *aodH = dynamic_cast<AliAODInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());\r
136         \r
137     if (!aodH) {\r
138       TObject *handler = AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler();\r
139       AliInfo("Has output handler ");\r
140       if( handler && handler->InheritsFrom("AliAODHandler") ) {\r
141         AliInfo("Selected AOD analysis");\r
142 \r
143         fAOD = ((AliAODHandler*)handler)->GetAOD();\r
144         fAnalysisType = 2;\r
145       }\r
146       else {\r
147         AliWarning("Selected AOD reader but no AOD handler found");\r
148       }\r
149     } \r
150     else {\r
151       AliInfo("Selected AOD analysis");\r
152       fAnalysisType = 2;\r
153       \r
154       fAOD = aodH->GetEvent();\r
155     }\r
156   }\r
157 \r
158   if ((!fAOD) && (!fESD)) {\r
159     AliWarning("Wrong analysis type: Only ESD and AOD types are allowed!");\r
160   }\r
161 }\r
162 \r
163 //________________________________________________________________________\r
164 void AliAnalysisTaskFemto::CreateOutputObjects() {\r
165   AliInfo("Creating Femto Analysis objects\n");\r
166 \r
167   gSystem->SetIncludePath("-I$ROOTSYS/include -I./STEERBase/ -I./ESD/ -I./AOD/ -I./ANALYSIS/ -I./ANALYSISalice/ -I./PWG2AOD/AOD -I./PWG2femtoscopy/FEMTOSCOPY/AliFemto -I./PWG2femtoscopyUser/FEMTOSCOPY/AliFemtoUser");\r
168   //  char fcm[2000];\r
169 //   sprintf(fcm, "%s++", fConfigMacro);\r
170 //   gROOT->LoadMacro(fcm);\r
171   gROOT->LoadMacro(fConfigMacro);\r
172   //  fJetFinder = (AliJetFinder*) gInterpreter->ProcessLine("ConfigJetAnalysis()");\r
173   SetFemtoManager((AliFemtoManager *) gInterpreter->ProcessLine("ConfigFemtoAnalysis()"));\r
174 \r
175   TList *tOL;\r
176   fOutputList = fManager->Analysis(0)->GetOutputList();\r
177 \r
178   for (unsigned int ian = 1; ian<fManager->AnalysisCollection()->size(); ian++) {\r
179     tOL = fManager->Analysis(ian)->GetOutputList();\r
180 \r
181     TIter nextListCf(tOL);\r
182     while (TObject *obj = nextListCf()) {\r
183       fOutputList->Add(obj);\r
184     }\r
185 \r
186     delete tOL;\r
187   }\r
188 }\r
189 \r
190 //________________________________________________________________________\r
191 void AliAnalysisTaskFemto::Exec(Option_t *) {\r
192   // Task making a femtoscopic analysis.\r
193 \r
194   if (fAnalysisType==1) {\r
195     if (!fESD) {\r
196       AliWarning("fESD not available");\r
197       return;\r
198     }\r
199 \r
200     //Get MC data\r
201     AliMCEventHandler*    mctruth = (AliMCEventHandler*) \r
202       ((AliAnalysisManager::GetAnalysisManager())->GetMCtruthEventHandler());\r
203     \r
204     AliGenHijingEventHeader *hdh = 0;\r
205     if(mctruth) {\r
206       fStack = mctruth->MCEvent()->Stack();\r
207 \r
208       AliGenCocktailEventHeader *hd = dynamic_cast<AliGenCocktailEventHeader *> (mctruth->MCEvent()->GenEventHeader());\r
209       \r
210       if (hd) {\r
211         \r
212         //      AliInfo ("Got MC cocktail event header %p\n", (void *) hd);\r
213         TList *lhd = hd->GetHeaders();\r
214         //      AliInfo ("Got list of headers %d\n", lhd->GetEntries());\r
215         \r
216         for (int iterh=0; iterh<lhd->GetEntries(); iterh++) \r
217           {\r
218             hdh = dynamic_cast<AliGenHijingEventHeader *> (lhd->At(iterh));\r
219             //      AliInfo ("HIJING header at %i is %p\n", iterh, (void *) hdh);\r
220           }\r
221       }    \r
222     }\r
223 \r
224     // Get ESD\r
225     AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());\r
226     \r
227     if (!esdH) {\r
228       AliWarning("Could not get ESDInputHandler");\r
229       return;\r
230     } \r
231     else {\r
232       fESD = esdH->GetEvent();\r
233     }\r
234 \r
235     AliInfo(Form("Tracks in ESD: %d \n",fESD->GetNumberOfTracks()));\r
236 \r
237     if (fESD->GetNumberOfTracks() >= 0) {\r
238     \r
239       if (!fReader) {\r
240         AliWarning("No ESD reader for ESD analysis !\n");\r
241       }\r
242       \r
243       AliFemtoEventReaderESDChain* fesdc = dynamic_cast<AliFemtoEventReaderESDChain *> (fReader);\r
244       if (fesdc)\r
245         {\r
246           // Process the event with no Kine information\r
247           fesdc->SetESDSource(fESD);\r
248           fManager->ProcessEvent();\r
249         }\r
250       AliFemtoEventReaderESDChainKine* fesdck = dynamic_cast<AliFemtoEventReaderESDChainKine *> (fReader);\r
251       if (fesdck) \r
252         {\r
253           // Process the event with Kine information\r
254           fesdck->SetESDSource(fESD);\r
255           fesdck->SetStackSource(fStack);\r
256           \r
257           fesdck->SetGenEventHeader(hdh);\r
258           fManager->ProcessEvent();\r
259         }\r
260       AliFemtoEventReaderStandard* fstd = dynamic_cast<AliFemtoEventReaderStandard *> (fReader);\r
261       if (fstd) \r
262         {\r
263           // Process the event with Kine information\r
264           fstd->SetESDSource(fESD);\r
265           if (mctruth) {\r
266             fstd->SetStackSource(fStack);\r
267             fstd->SetGenEventHeader(hdh);\r
268             fstd->SetInputType(AliFemtoEventReaderStandard::kESDKine);\r
269           }\r
270           else\r
271             fstd->SetInputType(AliFemtoEventReaderStandard::kESD);\r
272           fManager->ProcessEvent();\r
273         }\r
274     } \r
275 \r
276     // Post the output histogram list\r
277     PostData(0, fOutputList);\r
278   }\r
279 \r
280   if (fAnalysisType==2) {    \r
281     if (!fAOD) {\r
282       AliWarning("fAOD not available");\r
283       return;\r
284     }\r
285 \r
286     // Get AOD\r
287 //     AliAODInputHandler *aodH = dynamic_cast<AliAODInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());\r
288       \r
289 //     if (!aodH) {\r
290 //       AliWarning("Could not get AODInputHandler");\r
291 //       return;\r
292 //     } \r
293 //     else {\r
294 \r
295 //       fAOD = aodH->GetEvent();\r
296 //     }\r
297 \r
298     AliInfo(Form("Tracks in AOD: %d \n",fAOD->GetNumberOfTracks()));\r
299     \r
300     if (fAOD->GetNumberOfTracks() > 0) {\r
301       if (!fReader) {\r
302         AliWarning("No AOD reader for AOD analysis! \n");\r
303       }\r
304       else {\r
305         AliFemtoEventReaderAODChain* faodc = dynamic_cast<AliFemtoEventReaderAODChain *> (fReader);\r
306 \r
307         if (faodc) {\r
308           // Process the event\r
309           faodc->SetAODSource(fAOD);\r
310           fManager->ProcessEvent();\r
311         }\r
312         AliFemtoEventReaderStandard* fstd = dynamic_cast<AliFemtoEventReaderStandard *> (fReader);\r
313 \r
314         if (fstd) {\r
315           // Process the event\r
316           fstd->SetAODSource(fAOD);\r
317           fstd->SetInputType(AliFemtoEventReaderStandard::kAOD);\r
318           fManager->ProcessEvent();\r
319         }\r
320       }\r
321     } \r
322 \r
323     // Post the output histogram list\r
324     PostData(0, fOutputList);\r
325   }\r
326 }      \r
327 \r
328 //________________________________________________________________________\r
329 void AliAnalysisTaskFemto::Terminate(Option_t *) {\r
330   // Do the final processing\r
331   if (fManager) {\r
332     fManager->Finish();\r
333   }\r
334 }\r
335 //________________________________________________________________________\r
336 void AliAnalysisTaskFemto:: FinishTaskOutput() {\r
337   // Do the final processing\r
338   if (fManager) {\r
339     fManager->Finish();\r
340   }\r
341 }\r
342 //________________________________________________________________________\r
343 void AliAnalysisTaskFemto::SetFemtoReaderESD(AliFemtoEventReaderESDChain *aReader)\r
344 {\r
345   AliInfo("Selecting Femto reader for ESD\n");\r
346   fReader = aReader;\r
347 }\r
348 //________________________________________________________________________\r
349 void AliAnalysisTaskFemto::SetFemtoReaderESDKine(AliFemtoEventReaderESDChainKine *aReader)\r
350 {\r
351   AliInfo("Selecting Femto reader for ESD with Kinematics information\n");\r
352   fReader = aReader;\r
353 }\r
354 //________________________________________________________________________\r
355 void AliAnalysisTaskFemto::SetFemtoReaderAOD(AliFemtoEventReaderAODChain *aReader)\r
356 {\r
357   AliInfo("Selecting Femto reader for AOD\n");\r
358   fReader = aReader;\r
359 }\r
360 void AliAnalysisTaskFemto::SetFemtoReaderStandard(AliFemtoEventReaderStandard *aReader)\r
361 {\r
362   AliInfo("Selecting Standard all-purpose Femto reader\n");\r
363   fReader = aReader;\r
364 }\r
365 //________________________________________________________________________\r
366 void AliAnalysisTaskFemto::SetFemtoManager(AliFemtoManager *aManager)\r
367 {\r
368   fManager = aManager;\r
369   AliInfo(Form("Got reader %p\n", (void *) aManager->EventReader()));\r
370   AliFemtoEventReaderESDChain     *tReaderESDChain     = dynamic_cast<AliFemtoEventReaderESDChain *> (aManager->EventReader());\r
371   AliFemtoEventReaderESDChainKine *tReaderESDChainKine = dynamic_cast<AliFemtoEventReaderESDChainKine *> (aManager->EventReader());\r
372   AliFemtoEventReaderAODChain     *tReaderAODChain     = dynamic_cast<AliFemtoEventReaderAODChain *> (aManager->EventReader());\r
373   AliFemtoEventReaderStandard     *tReaderStandard     = dynamic_cast<AliFemtoEventReaderStandard *> (aManager->EventReader());\r
374 \r
375   if ((!tReaderESDChain) && (!tReaderESDChainKine) && (!tReaderAODChain) && (!tReaderStandard)) {\r
376     AliWarning("No AliFemto event reader created. Will not run femto analysis.\n");\r
377     return;\r
378   }\r
379   if (tReaderESDChain) SetFemtoReaderESD(tReaderESDChain);\r
380   if (tReaderESDChainKine) SetFemtoReaderESDKine(tReaderESDChainKine);\r
381   if (tReaderAODChain) SetFemtoReaderAOD(tReaderAODChain);\r
382   if (tReaderStandard) SetFemtoReaderStandard(tReaderStandard);\r
383 }\r
384 \r