]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FEMTOSCOPY/AliFemto/AliAnalysisTaskFemto.cxx
Select configuration via argumen
[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   printf("   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     Printf("ERROR: Could not read chain from input slot 0");\r
114   } \r
115   else {\r
116     AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());\r
117     \r
118     if(esdH) {\r
119       cout << "Selected ESD analysis" << endl;\r
120       fAnalysisType = 1;\r
121       \r
122       if (!esdH) {\r
123         Printf("ERROR: Could not get ESDInputHandler");\r
124       } \r
125       else {\r
126         fESD = esdH->GetEvent();\r
127       }\r
128     }\r
129     else {\r
130       AliAODInputHandler *aodH = dynamic_cast<AliAODInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());\r
131       \r
132       if (!aodH) {\r
133         Printf("ERROR: Could not get AODInputHandler");\r
134       } \r
135       else {\r
136         cout << "Selected AOD analysis" << endl;\r
137         fAnalysisType = 2;\r
138 \r
139         fAOD = aodH->GetEvent();\r
140       }\r
141     }\r
142     if ((!fAOD) && (!fESD)) {\r
143       Printf("Wrong analysis type: Only ESD and AOD types are allowed!");\r
144     }\r
145   }\r
146   \r
147   \r
148 }\r
149 \r
150 //________________________________________________________________________\r
151 void AliAnalysisTaskFemto::CreateOutputObjects() {\r
152   printf("Creating Femto Analysis objects\n");\r
153 \r
154   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
155   gROOT->LoadMacro(fConfigMacro);\r
156   //  fJetFinder = (AliJetFinder*) gInterpreter->ProcessLine("ConfigJetAnalysis()");\r
157   SetFemtoManager((AliFemtoManager *) gInterpreter->ProcessLine("ConfigFemtoAnalysis()"));\r
158 \r
159   TList *tOL;\r
160   fOutputList = fManager->Analysis(0)->GetOutputList();\r
161 \r
162   for (unsigned int ian = 1; ian<fManager->AnalysisCollection()->size(); ian++) {\r
163     tOL = fManager->Analysis(ian)->GetOutputList();\r
164 \r
165     TIter nextListCf(tOL);\r
166     while (TObject *obj = nextListCf()) {\r
167       fOutputList->Add(obj);\r
168     }\r
169 \r
170     delete tOL;\r
171   }\r
172 }\r
173 \r
174 //________________________________________________________________________\r
175 void AliAnalysisTaskFemto::Exec(Option_t *) {\r
176   // Task making a femtoscopic analysis.\r
177 \r
178   if (fAnalysisType==1) {\r
179     if (!fESD) {\r
180       Printf("ERROR: fESD not available");\r
181       return;\r
182     }\r
183 \r
184     //Get MC data\r
185     AliMCEventHandler*    mctruth = (AliMCEventHandler*) \r
186       ((AliAnalysisManager::GetAnalysisManager())->GetMCtruthEventHandler());\r
187     \r
188     AliGenHijingEventHeader *hdh = 0;\r
189     if(mctruth) {\r
190       fStack = mctruth->MCEvent()->Stack();\r
191 \r
192       AliGenCocktailEventHeader *hd = dynamic_cast<AliGenCocktailEventHeader *> (mctruth->MCEvent()->GenEventHeader());\r
193       \r
194       if (hd) {\r
195         \r
196         //      printf ("Got MC cocktail event header %p\n", (void *) hd);\r
197         TList *lhd = hd->GetHeaders();\r
198         //      printf ("Got list of headers %d\n", lhd->GetEntries());\r
199         \r
200         for (int iterh=0; iterh<lhd->GetEntries(); iterh++) \r
201           {\r
202             hdh = dynamic_cast<AliGenHijingEventHeader *> (lhd->At(iterh));\r
203             //      printf ("HIJING header at %i is %p\n", iterh, (void *) hdh);\r
204           }\r
205       }    \r
206     }\r
207 \r
208     // Get ESD\r
209     AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());\r
210     \r
211     if (!esdH) {\r
212       Printf("ERROR: Could not get ESDInputHandler");\r
213       return;\r
214     } \r
215     else {\r
216       fESD = esdH->GetEvent();\r
217     }\r
218 \r
219     printf("Tracks in ESD: %d \n",fESD->GetNumberOfTracks());\r
220 \r
221     if (fESD->GetNumberOfTracks() >= 0) {\r
222     \r
223       if (!fReader) {\r
224         printf("ERROR: No ESD reader for ESD analysis !\n");\r
225       }\r
226       \r
227       AliFemtoEventReaderESDChain* fesdc = dynamic_cast<AliFemtoEventReaderESDChain *> (fReader);\r
228       if (fesdc)\r
229         {\r
230           // Process the event with no Kine information\r
231           fesdc->SetESDSource(fESD);\r
232           fManager->ProcessEvent();\r
233         }\r
234       AliFemtoEventReaderESDChainKine* fesdck = dynamic_cast<AliFemtoEventReaderESDChainKine *> (fReader);\r
235       if (fesdck) \r
236         {\r
237           // Process the event with Kine information\r
238           fesdck->SetESDSource(fESD);\r
239           fesdck->SetStackSource(fStack);\r
240           \r
241           fesdck->SetGenEventHeader(hdh);\r
242           fManager->ProcessEvent();\r
243         }\r
244     } \r
245 \r
246     // Post the output histogram list\r
247     PostData(0, fOutputList);\r
248   }\r
249   \r
250   if (fAnalysisType==2) {    \r
251     if (!fAOD) {\r
252       Printf("ERROR: fAOD not available");\r
253       return;\r
254     }\r
255 \r
256     // Get AOD\r
257     AliAODInputHandler *aodH = dynamic_cast<AliAODInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());\r
258       \r
259     if (!aodH) {\r
260       Printf("ERROR: Could not get AODInputHandler");\r
261       return;\r
262     } \r
263     else {\r
264 \r
265       fAOD = aodH->GetEvent();\r
266     }\r
267 \r
268     printf("Tracks in AOD: %d \n",fAOD->GetNumberOfTracks());\r
269     \r
270     if (fAOD->GetNumberOfTracks() > 0) {\r
271       if (!fReader) {\r
272         printf("ERROR: No AOD reader for AOD analysis! \n");\r
273       }\r
274       else {\r
275         AliFemtoEventReaderAODChain* faodc = dynamic_cast<AliFemtoEventReaderAODChain *> (fReader);\r
276 \r
277         if (faodc) {\r
278           // Process the event\r
279           faodc->SetAODSource(fAOD);\r
280           fManager->ProcessEvent();\r
281         }\r
282       }\r
283     } \r
284 \r
285     // Post the output histogram list\r
286     PostData(0, fOutputList);\r
287   }\r
288 }      \r
289 \r
290 //________________________________________________________________________\r
291 void AliAnalysisTaskFemto::Terminate(Option_t *) {\r
292   // Do the final processing\r
293   if (fManager) {\r
294     fManager->Finish();\r
295   }\r
296 }\r
297 //________________________________________________________________________\r
298 void AliAnalysisTaskFemto:: FinishTaskOutput() {\r
299   // Do the final processing\r
300   if (fManager) {\r
301     fManager->Finish();\r
302   }\r
303 }\r
304 //________________________________________________________________________\r
305 void AliAnalysisTaskFemto::SetFemtoReaderESD(AliFemtoEventReaderESDChain *aReader)\r
306 {\r
307   printf("Selectring Femto reader for ESD\n");\r
308   fReader = aReader;\r
309 }\r
310 //________________________________________________________________________\r
311 void AliAnalysisTaskFemto::SetFemtoReaderESDKine(AliFemtoEventReaderESDChainKine *aReader)\r
312 {\r
313   printf("Selectring Femto reader for ESD with Kinematics information\n");\r
314   fReader = aReader;\r
315 }\r
316 //________________________________________________________________________\r
317 void AliAnalysisTaskFemto::SetFemtoReaderAOD(AliFemtoEventReaderAODChain *aReader)\r
318 {\r
319   printf("Selecting Femto reader for AOD\n");\r
320   fReader = aReader;\r
321 }\r
322 //________________________________________________________________________\r
323 void AliAnalysisTaskFemto::SetFemtoManager(AliFemtoManager *aManager)\r
324 {\r
325   fManager = aManager;\r
326   printf("Got reader %p\n", (void *) aManager->EventReader());\r
327   AliFemtoEventReaderESDChain     *tReaderESDChain     = dynamic_cast<AliFemtoEventReaderESDChain *> (aManager->EventReader());\r
328   AliFemtoEventReaderESDChainKine *tReaderESDChainKine = dynamic_cast<AliFemtoEventReaderESDChainKine *> (aManager->EventReader());\r
329   AliFemtoEventReaderAODChain     *tReaderAODChain     = dynamic_cast<AliFemtoEventReaderAODChain *> (aManager->EventReader());\r
330 \r
331   if ((!tReaderESDChain) && (!tReaderESDChainKine) && (!tReaderAODChain)) {\r
332     printf("No AliFemto event reader created. Will not run femto analysis.\n");\r
333     return;\r
334   }\r
335   if (tReaderESDChain) SetFemtoReaderESD(tReaderESDChain);\r
336   if (tReaderESDChainKine) SetFemtoReaderESDKine(tReaderESDChainKine);\r
337   if (tReaderAODChain) SetFemtoReaderAOD(tReaderAODChain);\r
338 }\r
339 \r