]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/CaloTrackCorrBase/AliAnaCaloTrackCorrMaker.cxx
AddTask added
[u/mrichter/AliRoot.git] / PWG / CaloTrackCorrBase / AliAnaCaloTrackCorrMaker.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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 // Steering class for particle (gamma, hadron) identification and correlation 
18 // analysis. It is called by the task class AliAnalysisTaskCaloTrackCorrelation 
19 // and it connects the input (ESD/AOD/MonteCarlo) got with AliCaloTrackReader 
20 // (produces TClonesArrays of AODs (TParticles in MC case if requested)), with 
21 // the analysis classes that derive from AliAnaCaloTrackCorrBaseClass
22 //
23 // -- Author: Gustavo Conesa (INFN-LNF, LPSC-Grenoble)
24
25 #include <cstdlib>
26
27 // --- ROOT system ---
28 #include "TClonesArray.h"
29 #include "TList.h"
30 #include "TH1I.h"
31 //#include <TObjectTable.h>
32
33 //---- AliRoot system ---- 
34 #include "AliAnalysisManager.h"
35 #include "AliInputEventHandler.h"
36 #include "AliAnaCaloTrackCorrBaseClass.h" 
37 #include "AliAnaCaloTrackCorrMaker.h" 
38
39 ClassImp(AliAnaCaloTrackCorrMaker)
40
41
42 //__________________________________________________
43 AliAnaCaloTrackCorrMaker::AliAnaCaloTrackCorrMaker() : 
44 TObject(),
45 fReader(0),                   fCaloUtils(0),
46 fOutputContainer(new TList ), fAnalysisContainer(new TList ),
47 fMakeHisto(kFALSE),           fMakeAOD(kFALSE), 
48 fAnaDebug(0),                 fCuts(new TList), 
49 fScaleFactor(-1),
50 fhNEvents(0),                 fhZVertex(0),                 
51 fhTrackMult(0),               fhCentrality(0),
52 fhEventPlaneAngle(0),
53 fhNMergedFiles(0),            fhScaleFactor(0)
54 {
55   //Default Ctor
56   if(fAnaDebug > 1 ) printf("*** Analysis Maker Constructor *** \n");
57   
58   //Initialize parameters, pointers and histograms
59   InitParameters();
60 }
61
62 //________________________________________________________________________________________
63 AliAnaCaloTrackCorrMaker::AliAnaCaloTrackCorrMaker(const AliAnaCaloTrackCorrMaker & maker) :   
64 TObject(),
65 fReader(),   //(new AliCaloTrackReader(*maker.fReader)),
66 fCaloUtils(),//(new AliCalorimeterUtils(*maker.fCaloUtils)),
67 fOutputContainer(new TList()), fAnalysisContainer(new TList()), 
68 fMakeHisto(maker.fMakeHisto),  fMakeAOD(maker.fMakeAOD),
69 fAnaDebug(maker.fAnaDebug),    fCuts(new TList()),
70 fScaleFactor(maker.fScaleFactor),
71 fhNEvents(maker.fhNEvents),    fhZVertex(maker.fhZVertex),    
72 fhTrackMult(maker.fhTrackMult),fhCentrality(maker.fhCentrality),
73 fhEventPlaneAngle(maker.fhEventPlaneAngle),
74 fhNMergedFiles(maker.fhNMergedFiles),          
75 fhScaleFactor(maker.fhScaleFactor)
76 {
77   // cpy ctor
78 }
79
80 //___________________________________________________
81 AliAnaCaloTrackCorrMaker::~AliAnaCaloTrackCorrMaker() 
82 {
83   // Remove all owned pointers.
84   
85   //  Do not delete it here, already done somewhere else, need to understand where.     
86   //  if (fOutputContainer) {
87   //    fOutputContainer->Clear();
88   //    delete fOutputContainer ;
89   //  }   
90   
91   if (fAnalysisContainer)
92   {
93     fAnalysisContainer->Delete();
94     delete fAnalysisContainer ;
95   }   
96   
97   if (fReader)    delete fReader ;
98   if (fCaloUtils) delete fCaloUtils ;
99   
100   if(fCuts)
101   {
102           fCuts->Delete();
103           delete fCuts;
104   }
105         
106 }
107
108 //__________________________________________________________________
109 void    AliAnaCaloTrackCorrMaker::AddAnalysis(TObject* ana, Int_t n) 
110 {
111   // Add analysis depending on AliAnaCaloTrackCorrBaseClass to list
112   
113   if ( fAnalysisContainer)
114   { 
115     fAnalysisContainer->AddAt(ana,n); 
116   }
117   else
118   { 
119     printf("AliAnaCaloTrackCorrMaker::AddAnalysis() - AnalysisContainer not initialized\n");
120     abort();
121   }
122 }  
123
124 //_________________________________________________________
125 TList * AliAnaCaloTrackCorrMaker::FillAndGetAODBranchList()
126
127         
128         // Get any new output AOD branches from analysis and put them in a list
129         // The list is filled in the maker, and new branch passed to the analysis frame
130         // AliAnalysisTaskCaloTrackCorrelation
131   
132         TList *aodBranchList = fReader->GetAODBranchList() ;
133   
134         for(Int_t iana = 0; iana <  fAnalysisContainer->GetEntries(); iana++)
135   {
136                 AliAnaCaloTrackCorrBaseClass * ana =  ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
137                 if(ana->NewOutputAOD()) aodBranchList->Add(ana->GetCreateOutputAODBranch());
138         }
139         
140         return aodBranchList ;
141         
142 }
143
144 //_______________________________________________________
145 TList * AliAnaCaloTrackCorrMaker::GetListOfAnalysisCuts()
146
147   
148         // Get the list of the cuts used for the analysis
149         // The list is filled in the maker, called by the task in LocalInit() and posted there
150   
151         for(Int_t iana = 0; iana <  fAnalysisContainer->GetEntries(); iana++)
152   {
153                 AliAnaCaloTrackCorrBaseClass * ana =  ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
154                 TObjString * objstring = ana->GetAnalysisCuts();
155     
156                 if(objstring)fCuts->Add(objstring);
157         }
158   
159         return fCuts ;
160   
161 }
162
163 //___________________________________________________
164 TList *AliAnaCaloTrackCorrMaker::GetOutputContainer()
165 {
166   // Fill the output list of histograms during the CreateOutputObjects stage.
167   
168   //Initialize calorimeters  geometry pointers
169   //GetCaloUtils()->InitPHOSGeometry();
170   //GetCaloUtils()->InitEMCALGeometry();
171   
172   //General event histograms
173   
174   fhNEvents      = new TH1I("hNEvents",   "Number of analyzed events"     , 1 , 0 , 1  ) ;
175   fhNEvents->SetYTitle("# events");
176   fOutputContainer->Add(fhNEvents);
177   
178   fhZVertex      = new TH1F("hZVertex", " Z vertex distribution"   , 200 , -50 , 50  ) ;
179   fhZVertex->SetXTitle("v_{z} (cm)");
180   fOutputContainer->Add(fhZVertex);
181   
182   fhTrackMult    = new TH1I("hTrackMult", "Number of tracks per events"   , 2000 , 0 , 2000  ) ;
183   fhTrackMult->SetXTitle("# tracks");
184   fOutputContainer->Add(fhTrackMult);
185   
186   fhCentrality   = new TH1F("hCentrality","Number of events in centrality bin",100,0.,100) ;
187   fhCentrality->SetXTitle("Centrality bin");
188   fOutputContainer->Add(fhCentrality) ;  
189   
190   fhEventPlaneAngle=new TH1F("hEventPlaneAngle","Number of events in event plane",100,0.,TMath::Pi()) ;
191   fhEventPlaneAngle->SetXTitle("EP angle (rad)");
192   fOutputContainer->Add(fhEventPlaneAngle) ;
193   
194   if(fScaleFactor > 0)
195   {
196     fhNMergedFiles = new TH1I("hNMergedFiles",   "Number of merged output files"     , 1 , 0 , 1  ) ;
197     fhNMergedFiles->SetYTitle("# files");
198     fhNMergedFiles->Fill(1); // Fill here with one entry, while merging it will count the rest
199     fOutputContainer->Add(fhNMergedFiles);
200     
201     fhScaleFactor = new TH1F("hScaleFactor",   "Number of merged output files"     , 1 , 0 , 1  ) ;
202     fhScaleFactor->SetYTitle("scale factor");
203     fhScaleFactor->SetBinContent(1,fScaleFactor); // Fill here 
204     fOutputContainer->Add(fhScaleFactor);    
205   }
206   
207   if(!fAnalysisContainer || fAnalysisContainer->GetEntries()==0)
208   {
209     printf("AliAnaCaloTrackCorrMaker::GetOutputContainer() - Analysis job list not initialized!!!\n");
210     return fOutputContainer;
211   }
212   
213   const Int_t buffersize = 255;
214   char newname[buffersize];
215   for(Int_t iana = 0; iana <  fAnalysisContainer->GetEntries(); iana++)
216   {
217     
218     AliAnaCaloTrackCorrBaseClass * ana =  ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
219     
220     if(fMakeHisto) // Analysis with histograms as output on
221     {
222       
223       //Fill container with appropriate histograms                      
224       TList * templist =  ana ->GetCreateOutputObjects(); 
225       templist->SetOwner(kFALSE); //Owner is fOutputContainer.
226       
227       for(Int_t i = 0; i < templist->GetEntries(); i++)
228       {
229         
230         //Add only  to the histogram name the name of the task
231         if(   strcmp((templist->At(i))->ClassName(),"TObjString")   ) 
232         {
233           snprintf(newname,buffersize, "%s%s", (ana->GetAddedHistogramsStringToName()).Data(), (templist->At(i))->GetName());  
234           ((TH1*) templist->At(i))->SetName(newname);
235         }
236         
237         //Add histogram to general container
238         fOutputContainer->Add(templist->At(i)) ;
239         
240       }
241       
242       delete templist;
243       
244     }// Analysis with histograms as output on
245     
246   }//Loop on analysis defined
247   
248   return fOutputContainer;
249   
250 }
251
252 //___________________________________
253 void AliAnaCaloTrackCorrMaker::Init()
254 {  
255   //Init container histograms and other common variables
256   // Fill the output list of histograms during the CreateOutputObjects stage.
257   
258   //Initialize reader
259   GetReader()->Init();
260   GetReader()->SetCaloUtils(GetCaloUtils()); // pass the calo utils pointer to the reader
261         
262   
263   if(!fAnalysisContainer || fAnalysisContainer->GetEntries()==0)
264   {
265     printf("AliAnaCaloTrackCorrMaker::GetOutputInit() - Analysis job list not initialized!!!\n");
266     return;
267   }
268   
269   for(Int_t iana = 0; iana <  fAnalysisContainer->GetEntries(); iana++)
270   {
271     
272     AliAnaCaloTrackCorrBaseClass * ana =  ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
273     
274     ana->SetReader(fReader);       // Set Reader for each analysis
275     ana->SetCaloUtils(fCaloUtils); // Set CaloUtils for each analysis
276     
277     ana->Init();
278     
279   }//Loop on analysis defined
280   
281 }
282
283 //_____________________________________________
284 void AliAnaCaloTrackCorrMaker::InitParameters()
285 {       
286   //Init data members
287   
288   fMakeHisto  = kTRUE;
289   fMakeAOD    = kTRUE; 
290   fAnaDebug   = 0; // No debugging info displayed by default
291         
292 }
293
294 //______________________________________________________________
295 void AliAnaCaloTrackCorrMaker::Print(const Option_t * opt) const
296 {       
297   //Print some relevant parameters set for the analysis
298         
299   if(! opt)
300     return;
301   
302   printf("***** Print: %s %s ******\n", GetName(), GetTitle() ) ;
303   printf("Debug level                =     %d\n", fAnaDebug   ) ;
304   printf("Produce Histo              =     %d\n", fMakeHisto  ) ;
305   printf("Produce AOD                =     %d\n", fMakeAOD    ) ;
306   printf("Number of analysis tasks   =     %d\n", fAnalysisContainer->GetEntries()) ;
307   
308   if(!strcmp("all",opt))
309   {
310           printf("Print analysis Tasks settings :\n") ;
311           for(Int_t iana = 0; iana<fAnalysisContainer->GetEntries(); iana++)
312     {
313                   ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana))->Print("");
314           }
315     
316           printf("Print analysis Reader settings :\n") ;
317           fReader->Print("");
318           printf("Print analysis Calorimeter Utils settings :\n") ;
319           fCaloUtils->Print("");
320     
321   }
322   
323
324
325 //_______________________________________________________________________
326 void AliAnaCaloTrackCorrMaker::ProcessEvent(const Int_t iEntry, 
327                                             const char * currentFileName)
328 {
329   //Process analysis for this event
330   
331   if(fMakeHisto && !fOutputContainer)
332   {
333     printf("AliAnaCaloTrackCorrMaker::ProcessEvent() - Histograms not initialized\n");
334     abort();
335   }
336         
337   if(fAnaDebug >= 0 )
338   { 
339                 printf("***  AliAnaCaloTrackCorrMaker::ProcessEvent() Event %d   ***  \n",iEntry);
340           if(fAnaDebug > 1 ) 
341     {
342                   printf("AliAnaCaloTrackCorrMaker::ProcessEvent() - Current File Name : %s\n", currentFileName);
343                   //printf("fAODBranchList %p, entries %d\n",fAODBranchList,fAODBranchList->GetEntries());
344           }
345   }
346   
347   //Each event needs an empty branch
348   TList * aodList = fReader->GetAODBranchList();
349   Int_t nAODBranches = aodList->GetEntries();
350   for(Int_t iaod = 0; iaod < nAODBranches; iaod++)
351   {
352           TClonesArray *tca = dynamic_cast<TClonesArray*> (aodList->At(iaod));
353           if(tca) tca->Clear("C");
354   }
355   
356   //Set geometry matrices before filling arrays, in case recalibration/position calculation etc is needed
357   fCaloUtils->AccessGeometry(fReader->GetInputEvent()); 
358   
359   //Set the AODB calibration, bad channels etc. parameters at least once
360   fCaloUtils->AccessOADB(fReader->GetInputEvent());     
361
362   
363   //Tell the reader to fill the data in the 3 detector lists
364   Bool_t ok = fReader->FillInputEvent(iEntry, currentFileName);
365   if(!ok)
366   {
367           if(fAnaDebug >= 1 )printf("*** Skip event *** %d \n",iEntry);
368           return ;
369   }
370         
371   //Magic line to write events to file
372   if(fReader->WriteDeltaAODToFile())AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler()->SetFillAOD(kTRUE);
373   
374   //printf(">>>>>>>>>> BEFORE >>>>>>>>>>>\n");
375   //gObjectTable->Print();
376   
377   //Loop on analysis algorithms
378   
379   if(fAnaDebug > 0 ) printf("*** Begin analysis *** \n");
380   
381   Int_t nana = fAnalysisContainer->GetEntries() ;
382   for(Int_t iana = 0; iana <  nana; iana++)
383   {
384     AliAnaCaloTrackCorrBaseClass * ana =  ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ; 
385     
386     ana->ConnectInputOutputAODBranches(); //Sets branches for each analysis
387     //Make analysis, create aods in aod branch or AODCaloClusters
388     if(fMakeAOD  )  ana->MakeAnalysisFillAOD()  ;
389     //Make further analysis with aod branch and fill histograms
390     if(fMakeHisto)  ana->MakeAnalysisFillHistograms()  ;
391     
392   }
393         
394   fReader->ResetLists();
395
396   // In case of mixing analysis, non triggered events are used,
397   // do not fill control histograms for a non requested triggered event
398   if(!fReader->IsEventTriggerAtSEOn())
399   {
400     AliAnalysisManager *manager = AliAnalysisManager::GetAnalysisManager();
401     AliInputEventHandler *inputHandler = dynamic_cast<AliInputEventHandler*>(manager->GetInputEventHandler());
402     
403     if(!inputHandler) return ;  
404     
405     UInt_t isTrigger = inputHandler->IsEventSelected() & fReader->GetEventTriggerMask();
406     if(!isTrigger) 
407     {
408       if(fAnaDebug > 0 ) printf("AliAnaCaloTrackMaker::ProcessEvent() - *** End analysis, MB for mixing *** \n");
409       
410       return;
411     }
412   }
413   
414   // Event control histograms
415  
416   fhNEvents        ->Fill(0); // Number of events analyzed
417   fhTrackMult      ->Fill(fReader->GetTrackMultiplicity()); 
418   fhCentrality     ->Fill(fReader->GetEventCentrality  ());
419   fhEventPlaneAngle->Fill(fReader->GetEventPlaneAngle  ());
420   
421   
422   Double_t v[3];
423   fReader->GetInputEvent()->GetPrimaryVertex()->GetXYZ(v) ;
424   fhZVertex->Fill(v[2]);
425   
426   
427   //printf(">>>>>>>>>> AFTER >>>>>>>>>>>\n");
428   //gObjectTable->Print();
429         
430   if(fAnaDebug > 0 ) printf("AliAnaCaloTrackMaker::ProcessEvent() - *** End analysis *** \n");
431   
432 }
433
434 //__________________________________________________________
435 void AliAnaCaloTrackCorrMaker::Terminate(TList * outputList)
436 {  
437   //Execute Terminate of analysis
438   //Do some final plots.
439   
440   if (!outputList) 
441   {
442           Error("Terminate", "No output list");
443           return;
444   }
445           
446   for(Int_t iana = 0; iana <  fAnalysisContainer->GetEntries(); iana++)
447   {
448     
449     AliAnaCaloTrackCorrBaseClass * ana =  ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
450     if(ana->MakePlotsOn())ana->Terminate(outputList);
451     
452   }//Loop on analysis defined
453   
454 }
455