]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/FEMTOSCOPY/AliFemto/AliAnalysisTaskFemto.cxx
Less bins in Mult monitors. New PID. Adapt to pPb analysis
[u/mrichter/AliRoot.git] / PWGCF / FEMTOSCOPY / AliFemto / AliAnalysisTaskFemto.cxx
1 //------------------------------------------------------
2 // AliAnalysisTaskFemto - A task for the analysis framework
3 // from the FEMTOSCOPY analysis of PWG2. Creates the necessary
4 // connection between the ESD or AOD input and the femtoscopic
5 // code.
6 // Author: Adam Kisiel, OSU; Adam.Kisiel@cern.ch
7 //------------------------------------------------------
8 #include "TROOT.h"
9 #include "TChain.h"
10 #include "TH1.h"
11 #include "TCanvas.h"
12 #include "TSystem.h"
13 #include "TFile.h"
14 #include "TInterpreter.h"
15
16 //#include "AliAnalysisTask.h"
17 #include "AliAnalysisTaskSE.h"
18
19 #include "AliESDEvent.h"
20
21 #include "AliFemtoAnalysis.h"
22 #include "AliAnalysisTaskFemto.h"
23 #include "AliVHeader.h"
24 #include "AliGenEventHeader.h"
25 #include "AliGenHijingEventHeader.h"
26 #include "AliGenCocktailEventHeader.h"
27
28 ClassImp(AliAnalysisTaskFemto)
29
30 // Default name for the setup macro of femto analysis  
31 // This function MUST be defined in the separate file !!!
32 // extern AliFemtoManager *ConfigFemtoAnalysis();
33
34 //________________________________________________________________________
35 AliAnalysisTaskFemto::AliAnalysisTaskFemto(const char *name, const char *aConfigMacro, const char *aConfigParams):
36 AliAnalysisTaskSE(name), //AliAnalysisTask(name,""), 
37   fESD(0), 
38   fESDpid(0),
39   fAOD(0),
40   fAODpidUtil(0),
41   fAODheader(0),
42   fStack(0),
43   fOutputList(0), 
44   fReader(0x0),
45   fManager(0x0),
46   fAnalysisType(0),
47     fConfigMacro(aConfigMacro),
48     fConfigParams(aConfigParams)
49 {
50   // Constructor.
51   // Input slot #0 works with an Ntuple
52   DefineInput(0, TChain::Class());
53   // Output slot #0 writes into a TH1 container
54   DefineOutput(0, TList::Class());
55
56 }
57 //________________________________________________________________________
58 AliAnalysisTaskFemto::AliAnalysisTaskFemto(const char *name, const char *aConfigMacro="ConfigFemtoAnalysis.C"): 
59     AliAnalysisTaskSE(name), //AliAnalysisTask(name,""), 
60     fESD(0), 
61     fESDpid(0),
62     fAOD(0),
63     fAODpidUtil(0),
64     fAODheader(0),
65     fStack(0),
66     fOutputList(0), 
67     fReader(0x0),
68     fManager(0x0),
69     fAnalysisType(0),
70     fConfigMacro(aConfigMacro),
71     fConfigParams("")
72 {
73   // Constructor.
74   // Input slot #0 works with an Ntuple
75   DefineInput(0, TChain::Class());
76   // Output slot #0 writes into a TH1 container
77   DefineOutput(0, TList::Class());
78
79 }
80
81 AliAnalysisTaskFemto::AliAnalysisTaskFemto(const AliAnalysisTaskFemto& aFemtoTask):
82     AliAnalysisTaskSE(aFemtoTask), //AliAnalysisTask(aFemtoTask), 
83     fESD(0), 
84     fESDpid(0),
85     fAOD(0),
86     fAODpidUtil(0),
87     fAODheader(0),
88     fStack(0),
89     fOutputList(0), 
90     fReader(0x0),
91     fManager(0x0),
92     fAnalysisType(0),
93     fConfigMacro(0),
94     fConfigParams(0)
95 {
96   // copy constructor
97   fESD = aFemtoTask.fESD; 
98   fESDpid = aFemtoTask.fESDpid; 
99   fAOD = aFemtoTask.fAOD; 
100   fAODpidUtil = aFemtoTask.fAODpidUtil;
101   fAODheader = aFemtoTask.fAODheader;
102   fStack = aFemtoTask.fStack;
103   fOutputList = aFemtoTask.fOutputList;   
104   fReader = aFemtoTask.fReader;       
105   fManager = aFemtoTask.fManager;      
106   fAnalysisType = aFemtoTask.fAnalysisType; 
107
108   fConfigMacro = aFemtoTask.fConfigMacro;
109   fConfigParams = aFemtoTask.fConfigParams;
110
111 }
112
113
114 AliAnalysisTaskFemto& AliAnalysisTaskFemto::operator=(const AliAnalysisTaskFemto& aFemtoTask){
115   // assignment operator
116   if (this == &aFemtoTask)
117     return *this;
118
119   fESD = aFemtoTask.fESD; 
120   fESDpid = aFemtoTask.fESDpid;
121   fAOD = aFemtoTask.fAOD; 
122   fAODpidUtil = aFemtoTask.fAODpidUtil;
123   fAODheader = aFemtoTask.fAODheader;
124   fStack = aFemtoTask.fStack;
125   fOutputList = aFemtoTask.fOutputList;   
126   fReader = aFemtoTask.fReader;       
127   fManager = aFemtoTask.fManager;      
128   fAnalysisType = aFemtoTask.fAnalysisType; 
129
130   fConfigMacro = aFemtoTask.fConfigMacro;
131   fConfigParams = aFemtoTask.fConfigParams;
132
133   return *this;
134 }
135
136 AliAnalysisTaskFemto::~AliAnalysisTaskFemto() 
137 {
138 }
139
140
141 //________________________________________________________________________
142 void AliAnalysisTaskFemto::ConnectInputData(Option_t *) {
143   AliInfo(Form("   ConnectInputData %s\n", GetName()));
144
145   fESD = 0;
146   fESDpid = 0;
147   fAOD = 0;
148   fAODpidUtil = 0;
149   fAODheader=0;
150   fAnalysisType = 0;
151
152   TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
153   if (!tree) {
154     AliWarning("Could not read chain from input slot 0");
155     return;
156   } 
157
158   AliFemtoEventReaderESDChain *femtoReader = dynamic_cast<AliFemtoEventReaderESDChain *> (fReader);
159   if ((dynamic_cast<AliFemtoEventReaderESDChain *> (fReader))) {
160     AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
161     
162     if(esdH) {
163       AliInfo("Selected ESD analysis");
164       fAnalysisType = 1;
165       
166 //       if (!esdH) {
167 //      AliWarning("Could not get ESDInputHandler");
168 //       } 
169 //       else {
170         fESD = esdH->GetEvent();
171         fESDpid = esdH->GetESDpid();
172         femtoReader->SetESDPid(fESDpid);
173 //       }
174     }
175 }
176   else if ((dynamic_cast<AliFemtoEventReaderKinematicsChain *> (fReader))) {
177     AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
178     
179     if(esdH) {
180       AliInfo("Selected ESD analysis");
181       fAnalysisType = 1;
182       
183 //       if (!esdH) {
184 //      AliWarning("Could not get ESDInputHandler");
185 //       } 
186 //       else {
187         fESD = esdH->GetEvent();
188         //fESDpid = esdH->GetESDpid();
189         //femtoReader->SetESDPid(fESDpid);
190 //       }
191     }
192   }
193
194
195   AliFemtoEventReaderESDChainKine *femtoReaderESDKine = dynamic_cast<AliFemtoEventReaderESDChainKine *> (fReader);
196   if ((dynamic_cast<AliFemtoEventReaderESDChainKine *> (fReader))) {
197     AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
198     
199     if(esdH) {
200       AliInfo("Selected ESD analysis");
201       fAnalysisType = 1;
202       
203 //       if (!esdH) {
204 //      AliWarning("Could not get ESDInputHandler");
205 //       } 
206 //       else {
207         fESD = esdH->GetEvent();
208         fESDpid = esdH->GetESDpid();
209         femtoReaderESDKine->SetESDPid(fESDpid);
210 //       }
211     }
212 }
213
214   //  AliFemtoEventReaderKinematicsChain *femtoReaderKine = dynamic_cast<AliFemtoEventReaderKinematicsChain *> (fReader);
215   if ((dynamic_cast<AliFemtoEventReaderKinematicsChain *> (fReader))) {
216     AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
217     
218     if(esdH) {
219       AliInfo("Selected ESD analysis");
220       fAnalysisType = 1;
221       
222       //       if (!esdH) {
223       //        AliWarning("Could not get ESDInputHandler");
224       //       } 
225       //       else {
226       fESD = esdH->GetEvent();
227       //fESDpid = esdH->GetESDpid();
228       //femtoReader->SetESDPid(fESDpid);
229       //       }
230     }
231   }
232   
233   
234   AliFemtoEventReaderAODChain *femtoReaderAOD = dynamic_cast<AliFemtoEventReaderAODChain *> (fReader);
235   if (dynamic_cast<AliFemtoEventReaderAODChain *> (fReader)) {
236     AliAODInputHandler *aodH = dynamic_cast<AliAODInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
237     
238     if (!aodH) {
239       TObject *handler = AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler();
240       AliInfo("Has output handler ");
241       if( handler && handler->InheritsFrom("AliAODHandler") ) {
242         AliInfo("Selected AOD analysis");
243         
244         fAOD = ((AliAODHandler*)handler)->GetAOD();
245         fAnalysisType = 2;
246       }
247       else {
248         AliWarning("Selected AOD reader but no AOD handler found");
249       }
250     } 
251     else {
252       AliInfo("Selected AOD analysis");
253       fAnalysisType = 2;
254       
255       fAOD = aodH->GetEvent();
256
257       fAODpidUtil = aodH->GetAODpidUtil(); //correct way
258       //fAODpidUtil = new AliAODpidUtil(); //not correct way
259       //      printf("aodH->GetAODpidUtil(): %x",aodH->GetAODpidUtil());
260       cout<<"AliAnalysisTaskFemto::AodpidUtil:"<<fAODpidUtil<<endl;
261       femtoReaderAOD->SetAODpidUtil(fAODpidUtil);
262
263       fAODheader = fAOD->GetHeader();
264       femtoReaderAOD->SetAODheader(fAODheader);
265    
266     }
267   }
268
269   if ((!fAOD) && (!fESD)) {
270     AliWarning("Wrong analysis type: Only ESD and AOD types are allowed!");
271   }
272 }
273
274 //________________________________________________________________________
275 void AliAnalysisTaskFemto::CreateOutputObjects() {
276   AliInfo("Creating Femto Analysis objects\n");
277
278   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");
279   //  char fcm[2000];
280 //   sprintf(fcm, "%s++", fConfigMacro);
281 //   gROOT->LoadMacro(fcm);
282   gROOT->LoadMacro(fConfigMacro);
283   //  fJetFinder = (AliJetFinder*) gInterpreter->ProcessLine("ConfigJetAnalysis()");
284   if (!fConfigParams)
285     SetFemtoManager((AliFemtoManager *) gInterpreter->ProcessLine("ConfigFemtoAnalysis()"));
286   else
287       SetFemtoManager((AliFemtoManager *) gInterpreter->ProcessLine(Form("ConfigFemtoAnalysis(%s)", fConfigParams.Data())));
288
289   TList *tOL;
290   fOutputList = fManager->Analysis(0)->GetOutputList();
291
292   for (unsigned int ian = 1; ian<fManager->AnalysisCollection()->size(); ian++) {
293     tOL = fManager->Analysis(ian)->GetOutputList();
294
295     TIter nextListCf(tOL);
296     while (TObject *obj = nextListCf()) {
297       fOutputList->Add(obj);
298     }
299
300     delete tOL;
301   }
302
303   PostData(0, fOutputList);
304 }
305
306 //________________________________________________________________________
307 void AliAnalysisTaskFemto::Exec(Option_t *) {
308   // Task making a femtoscopic analysis.
309
310   if (fAnalysisType==1) {
311     if (!fESD) {
312       AliWarning("fESD not available");
313       return;
314     }
315
316     //Get MC data
317     AliMCEventHandler*    mctruth = (AliMCEventHandler*) 
318       ((AliAnalysisManager::GetAnalysisManager())->GetMCtruthEventHandler());
319     
320     AliGenHijingEventHeader *hdh = 0;
321     if(mctruth) {
322       fStack = mctruth->MCEvent()->Stack();
323
324       AliGenCocktailEventHeader *hd = dynamic_cast<AliGenCocktailEventHeader *> (mctruth->MCEvent()->GenEventHeader());
325       
326       if (hd) {
327         
328         //      AliInfo ("Got MC cocktail event header %p\n", (void *) hd);
329         TList *lhd = hd->GetHeaders();
330         //      AliInfo ("Got list of headers %d\n", lhd->GetEntries());
331         
332         for (int iterh=0; iterh<lhd->GetEntries(); iterh++) 
333           {
334             hdh = dynamic_cast<AliGenHijingEventHeader *> (lhd->At(iterh));
335             //      AliInfo ("HIJING header at %i is %p\n", iterh, (void *) hdh);
336           }
337       }    
338     }
339
340     // Get ESD
341     AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
342     
343     if (!esdH) {
344       AliWarning("Could not get ESDInputHandler");
345       return;
346     } 
347     else {
348       fESD = esdH->GetEvent();
349       fESDpid = esdH->GetESDpid();   
350     }
351
352     AliInfo(Form("Tracks in ESD: %d \n",fESD->GetNumberOfTracks()));
353
354     if (fESD->GetNumberOfTracks() >= 0) {
355     
356       if (!fReader) {
357         AliWarning("No ESD reader for ESD analysis !\n");
358       }
359       
360       AliFemtoEventReaderESDChain* fesdc = dynamic_cast<AliFemtoEventReaderESDChain *> (fReader);
361       if (fesdc)
362         {
363           // Process the event with no Kine information
364           fesdc->SetESDSource(fESD);
365           fManager->ProcessEvent();
366         }
367
368         AliFemtoEventReaderKinematicsChain* fkinec = dynamic_cast<AliFemtoEventReaderKinematicsChain *> (fReader);
369       if (fkinec)
370         {
371           // Process the event with Kine information only
372           fkinec->SetStackSource(fStack);
373           fManager->ProcessEvent();
374         }
375
376
377       AliFemtoEventReaderESDChainKine* fesdck = dynamic_cast<AliFemtoEventReaderESDChainKine *> (fReader);
378       if (fesdck) 
379         {
380           // Process the event with Kine information
381           fesdck->SetESDSource(fESD);
382           fesdck->SetStackSource(fStack);
383           
384           fesdck->SetGenEventHeader(hdh);
385           fManager->ProcessEvent();
386         }
387       AliFemtoEventReaderStandard* fstd = dynamic_cast<AliFemtoEventReaderStandard *> (fReader);
388       if (fstd) 
389         {
390           // Process the event with Kine information
391           fstd->SetESDSource(fESD);
392           if (mctruth) {
393             fstd->SetStackSource(fStack);
394             fstd->SetGenEventHeader(hdh);
395             fstd->SetInputType(AliFemtoEventReaderStandard::kESDKine);
396           }
397           else
398             fstd->SetInputType(AliFemtoEventReaderStandard::kESD);
399           fManager->ProcessEvent();
400         }
401     } 
402
403     // Post the output histogram list
404     PostData(0, fOutputList);
405   }
406
407   if (fAnalysisType==2) {    
408     if (!fAOD) {
409       AliWarning("fAOD not available");
410       return;
411     }
412
413     // Get AOD
414 //     AliAODInputHandler *aodH = dynamic_cast<AliAODInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
415       
416 //     if (!aodH) {
417 //       AliWarning("Could not get AODInputHandler");
418 //       return;
419 //     } 
420 //     else {
421
422 //       fAOD = aodH->GetEvent();
423 //     }
424
425
426
427     Bool_t isSelected = (((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected() & fOfflineTriggerMask);
428   if(!isSelected) {//cout << "AliAnalysisTaskFemto: is not selected" << endl; 
429     return;}
430
431
432
433
434      AliInfo(Form("Tracks in AOD: %d \n",fAOD->GetNumberOfTracks()));
435     
436     if (fAOD->GetNumberOfTracks() > 0) {
437       if (!fReader) {
438         AliWarning("No AOD reader for AOD analysis! \n");
439       }
440       else {
441         AliFemtoEventReaderAODChain* faodc = dynamic_cast<AliFemtoEventReaderAODChain *> (fReader);
442
443         if (faodc) {
444           // Process the event
445           faodc->SetAODSource(fAOD);
446           fManager->ProcessEvent();
447         }
448         AliFemtoEventReaderStandard* fstd = dynamic_cast<AliFemtoEventReaderStandard *> (fReader);
449
450         if (fstd) {
451           // Process the event
452           fstd->SetAODSource(fAOD);
453           fstd->SetInputType(AliFemtoEventReaderStandard::kAOD);
454           fManager->ProcessEvent();
455         }
456       }
457     } 
458
459     // Post the output histogram list
460     PostData(0, fOutputList);
461   }
462 }      
463
464 //________________________________________________________________________
465 void AliAnalysisTaskFemto::Terminate(Option_t *) {
466   // Do the final processing
467   if (fManager) {
468     fManager->Finish();
469   }
470 }
471 //________________________________________________________________________
472 void AliAnalysisTaskFemto:: FinishTaskOutput() {
473   // Do the final processing
474   if (fManager) {
475     fManager->Finish();
476   }
477 }
478 //________________________________________________________________________
479 void AliAnalysisTaskFemto::SetFemtoReaderESD(AliFemtoEventReaderESDChain *aReader)
480 {
481   AliInfo("Selecting Femto reader for ESD\n");
482   fReader = aReader;
483 }
484 //________________________________________________________________________
485 void AliAnalysisTaskFemto::SetFemtoReaderESDKine(AliFemtoEventReaderESDChainKine *aReader)
486 {
487   AliInfo("Selecting Femto reader for ESD with Kinematics information\n");
488   fReader = aReader;
489 }
490 //________________________________________________________________________
491 void AliAnalysisTaskFemto::SetFemtoReaderAOD(AliFemtoEventReaderAODChain *aReader)
492 {
493   AliInfo("Selecting Femto reader for AOD\n");
494   fReader = aReader;
495 }
496 void AliAnalysisTaskFemto::SetFemtoReaderStandard(AliFemtoEventReaderStandard *aReader)
497 {
498   AliInfo("Selecting Standard all-purpose Femto reader\n");
499   fReader = aReader;
500 }
501 void AliAnalysisTaskFemto::SetFemtoReaderKinematics(AliFemtoEventReaderKinematicsChain *aReader)
502 {
503   printf("Selecting Femto reader for Kinematics (Monte Carlo) information\n");
504   fReader = aReader;
505 }
506 //________________________________________________________________________
507 void AliAnalysisTaskFemto::SetFemtoManager(AliFemtoManager *aManager)
508 {
509   fManager = aManager;
510   AliInfo(Form("Got reader %p\n", (void *) aManager->EventReader()));
511   AliFemtoEventReaderESDChain     *tReaderESDChain     = dynamic_cast<AliFemtoEventReaderESDChain *> (aManager->EventReader());
512   AliFemtoEventReaderESDChainKine *tReaderESDChainKine = dynamic_cast<AliFemtoEventReaderESDChainKine *> (aManager->EventReader());
513   AliFemtoEventReaderAODChain     *tReaderAODChain     = dynamic_cast<AliFemtoEventReaderAODChain *> (aManager->EventReader());
514   AliFemtoEventReaderStandard     *tReaderStandard     = dynamic_cast<AliFemtoEventReaderStandard *> (aManager->EventReader());
515   AliFemtoEventReaderKinematicsChain *tReaderKineChain = dynamic_cast<AliFemtoEventReaderKinematicsChain *> (aManager->EventReader());
516
517  if ((!tReaderESDChain) && (!tReaderESDChainKine) && (!tReaderAODChain) && (!tReaderStandard) && (!tReaderKineChain)) {
518     AliWarning("No AliFemto event reader created. Will not run femto analysis.\n");
519     return;
520   }
521   if (tReaderESDChain) SetFemtoReaderESD(tReaderESDChain);
522   if (tReaderESDChainKine) SetFemtoReaderESDKine(tReaderESDChainKine);
523   if (tReaderAODChain) SetFemtoReaderAOD(tReaderAODChain);
524   if (tReaderStandard) SetFemtoReaderStandard(tReaderStandard);
525   if (tReaderKineChain) SetFemtoReaderKinematics(tReaderKineChain);
526 }
527