]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/CaloTrackCorrBase/AliAnaCaloTrackCorrMaker.cxx
Coverity fixes (Diego)
[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 "Riostream.h"
32 //#include <TObjectTable.h>
33
34 //---- AliRoot system ---- 
35 #include "AliAnalysisManager.h"
36 #include "AliVEventHandler.h"
37 #include "AliAnaCaloTrackCorrBaseClass.h" 
38 #include "AliAnaCaloTrackCorrMaker.h" 
39
40 ClassImp(AliAnaCaloTrackCorrMaker)
41
42
43 //__________________________________________________
44 AliAnaCaloTrackCorrMaker::AliAnaCaloTrackCorrMaker() : 
45 TObject(),
46 fOutputContainer(new TList ), fAnalysisContainer(new TList ),
47 fMakeHisto(kFALSE),           fMakeAOD(kFALSE), 
48 fAnaDebug(0),                 fReader(0), 
49 fCaloUtils(0),                fCuts(new TList), 
50 fhNEvents(0x0),               fhTrackMult(0x0)
51 {
52   //Default Ctor
53   if(fAnaDebug > 1 ) printf("*** Analysis Maker Constructor *** \n");
54   
55   //Initialize parameters, pointers and histograms
56   InitParameters();
57 }
58
59 //________________________________________________________________________________________
60 AliAnaCaloTrackCorrMaker::AliAnaCaloTrackCorrMaker(const AliAnaCaloTrackCorrMaker & maker) :   
61 TObject(),
62 fOutputContainer(new TList()), fAnalysisContainer(new TList()), 
63 fMakeHisto(maker.fMakeHisto),  fMakeAOD(maker.fMakeAOD),
64 fAnaDebug(maker.fAnaDebug),    fReader(),//new AliCaloTrackReader(*maker.fReader)), 
65 fCaloUtils(),//(new AliCalorimeterUtils(*maker.fCaloUtils)),
66 fCuts(new TList()),
67 fhNEvents(maker.fhNEvents),    fhTrackMult(maker.fhTrackMult)
68 {
69   // cpy ctor
70 }
71
72 //___________________________________________________
73 AliAnaCaloTrackCorrMaker::~AliAnaCaloTrackCorrMaker() 
74 {
75   // Remove all owned pointers.
76   
77   //  Do not delete it here, already done somewhere else, need to understand where.     
78   //  if (fOutputContainer) {
79   //    fOutputContainer->Clear();
80   //    delete fOutputContainer ;
81   //  }   
82   
83   if (fAnalysisContainer) {
84     fAnalysisContainer->Delete();
85     delete fAnalysisContainer ;
86   }   
87   
88   if (fReader)    delete fReader ;
89   if (fCaloUtils) delete fCaloUtils ;
90   
91   if(fCuts){
92           fCuts->Delete();
93           delete fCuts;
94   }
95         
96 }
97
98 //_______________________________________________________
99 void    AliAnaCaloTrackCorrMaker::AddAnalysis(TObject* ana, Int_t n) 
100 {
101   // Add analysis depending on AliAnaCaloTrackCorrBaseClass to list
102   
103   if ( fAnalysisContainer)
104   { 
105     fAnalysisContainer->AddAt(ana,n); 
106   }
107   else
108   { 
109     printf("AliAnaCaloTrackCorrMaker::AddAnalysis() - AnalysisContainer not initialized\n");
110     abort();
111   }
112 }  
113
114 //_________________________________________________________
115 TList * AliAnaCaloTrackCorrMaker::FillAndGetAODBranchList()
116
117         
118         // Get any new output AOD branches from analysis and put them in a list
119         // The list is filled in the maker, and new branch passed to the analysis frame
120         // AliAnalysisTaskCaloTrackCorrelation
121   
122         TList *aodBranchList = fReader->GetAODBranchList() ;
123   
124         for(Int_t iana = 0; iana <  fAnalysisContainer->GetEntries(); iana++){
125                 
126                 AliAnaCaloTrackCorrBaseClass * ana =  ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
127                 if(ana->NewOutputAOD()) aodBranchList->Add(ana->GetCreateOutputAODBranch());
128     
129         }
130         
131         return aodBranchList ;
132         
133 }
134
135 //_______________________________________________________
136 TList * AliAnaCaloTrackCorrMaker::GetListOfAnalysisCuts()
137
138   
139         // Get the list of the cuts used for the analysis
140         // The list is filled in the maker, called by the task in LocalInit() and posted there
141
142   
143         for(Int_t iana = 0; iana <  fAnalysisContainer->GetEntries(); iana++){
144                 AliAnaCaloTrackCorrBaseClass * ana =  ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
145                 TObjString * objstring = ana->GetAnalysisCuts();
146
147                 if(objstring)fCuts->Add(objstring);
148         }
149
150         return fCuts ;
151   
152 }
153
154 //___________________________________________________
155 TList *AliAnaCaloTrackCorrMaker::GetOutputContainer()
156 {
157   // Fill the output list of histograms during the CreateOutputObjects stage.
158   
159   //Initialize calorimeters  geometry pointers
160   GetCaloUtils()->InitPHOSGeometry();
161   GetCaloUtils()->InitEMCALGeometry();
162   
163   if(!fAnalysisContainer || fAnalysisContainer->GetEntries()==0){
164     printf("AliAnaCaloTrackCorrMaker::GetOutputContainer() - Analysis job list not initialized!!!\n");
165     //abort();
166   }
167   else{
168     const Int_t buffersize = 255;
169     char newname[255];
170     for(Int_t iana = 0; iana <  fAnalysisContainer->GetEntries(); iana++){
171       AliAnaCaloTrackCorrBaseClass * ana =  ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
172       if(fMakeHisto){// Analysis with histograms as output on
173         //Fill container with appropriate histograms                    
174         TList * templist =  ana ->GetCreateOutputObjects(); 
175         templist->SetOwner(kFALSE); //Owner is fOutputContainer.
176         for(Int_t i = 0; i < templist->GetEntries(); i++){
177           
178           //Add only  to the histogram name the name of the task
179           if(   strcmp((templist->At(i))->ClassName(),"TObjString")   ) {
180             snprintf(newname,buffersize, "%s%s", (ana->GetAddedHistogramsStringToName()).Data(), (templist->At(i))->GetName());  
181             ((TH1*) templist->At(i))->SetName(newname);
182           }
183           //Add histogram to general container
184           fOutputContainer->Add(templist->At(i)) ;
185         }
186         delete templist;
187       }// Analysis with histograms as output on
188     }//Loop on analysis defined
189   }//Analysis list available
190   
191   //General event histograms
192   
193   fhNEvents        = new TH1I("hNEvents",   "Number of analyzed events"     , 1 , 0 , 1  ) ;
194   fOutputContainer->Add(fhNEvents);
195   fhTrackMult      = new TH1I("hTrackMult", "Number of tracks per events"   , 2000 , 0 , 2000  ) ;
196   fOutputContainer->Add(fhTrackMult);
197   
198   return fOutputContainer;
199   
200 }
201
202 //___________________________________
203 void AliAnaCaloTrackCorrMaker::Init()
204 {  
205   //Init container histograms and other common variables
206   // Fill the output list of histograms during the CreateOutputObjects stage.
207   
208   //Initialize reader
209   GetReader()->Init();
210   GetReader()->SetCaloUtils(GetCaloUtils()); // pass the calo utils pointer to the reader
211         
212   
213   if(!fAnalysisContainer || fAnalysisContainer->GetEntries()==0){
214     printf("AliAnaCaloTrackCorrMaker::GetOutputInit() - Analysis job list not initialized!!!\n");
215     //abort();
216   }
217         else{
218     
219     for(Int_t iana = 0; iana <  fAnalysisContainer->GetEntries(); iana++){
220       
221       AliAnaCaloTrackCorrBaseClass * ana =  ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
222       ana->SetReader(fReader); //SetReader for each analysis
223       ana->SetCaloUtils(fCaloUtils); //Set CaloUtils for each analysis
224       
225       ana->Init();
226       
227     }//Loop on analysis defined
228   }//Analysis list available
229 }
230
231 //_____________________________________________
232 void AliAnaCaloTrackCorrMaker::InitParameters()
233 {       
234   //Init data members
235   
236   fMakeHisto  = kTRUE;
237   fMakeAOD    = kTRUE; 
238   fAnaDebug   = 0; // No debugging info displayed by default
239         
240 }
241
242 //______________________________________________________________
243 void AliAnaCaloTrackCorrMaker::Print(const Option_t * opt) const
244 {       
245   //Print some relevant parameters set for the analysis
246         
247   if(! opt)
248     return;
249   
250   printf("***** Print: %s %s ******\n", GetName(), GetTitle() ) ;
251   printf("Debug level                =     %d\n", fAnaDebug) ;
252   printf("Produce Histo              =     %d\n", fMakeHisto) ;
253   printf("Produce AOD                =     %d\n", fMakeAOD) ;
254   printf("Number of analysis tasks   =     %d\n", fAnalysisContainer->GetEntries()) ;
255   
256   if(!strcmp("all",opt)){
257           printf("Print analysis Tasks settings :\n") ;
258           for(Int_t iana = 0; iana<fAnalysisContainer->GetEntries(); iana++){
259                   ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana))->Print("");
260           }
261     
262           printf("Print analysis Reader settings :\n") ;
263           fReader->Print("");
264           printf("Print analysis Calorimeter Utils settings :\n") ;
265           fCaloUtils->Print("");
266   }
267
268
269 //_______________________________________________________________________
270 void AliAnaCaloTrackCorrMaker::ProcessEvent(const Int_t iEntry, 
271                                             const char * currentFileName)
272 {
273   //Process analysis for this event
274   
275   if(fMakeHisto && !fOutputContainer){
276     printf("AliAnaCaloTrackCorrMaker::ProcessEvent() - Histograms not initialized\n");
277     abort();
278   }
279         
280   if(fAnaDebug >= 0 ){ 
281                 printf("***  AliAnaCaloTrackCorrMaker::ProcessEvent() Event %d   ***  \n",iEntry);
282           if(fAnaDebug > 1 ) {
283                   printf("AliAnaCaloTrackCorrMaker::ProcessEvent() - Current File Name : %s\n", currentFileName);
284                   //printf("fAODBranchList %p, entries %d\n",fAODBranchList,fAODBranchList->GetEntries());
285           }
286   }
287   
288   //Each event needs an empty branch
289   TList * aodList = fReader->GetAODBranchList();
290   Int_t nAODBranches = aodList->GetEntries();
291   for(Int_t iaod = 0; iaod < nAODBranches; iaod++){
292           TClonesArray *tca = dynamic_cast<TClonesArray*> (aodList->At(iaod));
293           if(tca) tca->Clear("C");
294   }
295   
296   //Set geometry matrices before filling arrays, in case recalibration/position calculation etc is needed
297   fCaloUtils->SetGeometryTransformationMatrices(fReader->GetInputEvent());      
298   
299   //Tell the reader to fill the data in the 3 detector lists
300   Bool_t ok = fReader->FillInputEvent(iEntry, currentFileName);
301   if(!ok){
302           if(fAnaDebug >= 1 )printf("*** Skip event *** %d \n",iEntry);
303           return ;
304   }
305         
306   //Magic line to write events to file
307   if(fReader->WriteDeltaAODToFile())AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler()->SetFillAOD(kTRUE);
308   
309   //printf(">>>>>>>>>> BEFORE >>>>>>>>>>>\n");
310   //gObjectTable->Print();
311   //Loop on analysis algorithms
312   if(fAnaDebug > 0 ) printf("*** Begin analysis *** \n");
313   Int_t nana = fAnalysisContainer->GetEntries() ;
314   for(Int_t iana = 0; iana <  nana; iana++){
315     AliAnaCaloTrackCorrBaseClass * ana =  ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ; 
316     
317     ana->ConnectInputOutputAODBranches(); //Sets branches for each analysis
318     //Make analysis, create aods in aod branch or AODCaloClusters
319     if(fMakeAOD  )  ana->MakeAnalysisFillAOD()  ;
320     //Make further analysis with aod branch and fill histograms
321     if(fMakeHisto)  ana->MakeAnalysisFillHistograms()  ;
322     
323   }
324         
325   fhNEvents->Fill(0); //Event analyzed
326   fhTrackMult->Fill(fReader->GetTrackMultiplicity()); 
327   
328   fReader->ResetLists();
329   
330   //printf(">>>>>>>>>> AFTER >>>>>>>>>>>\n");
331   //gObjectTable->Print();
332         
333   if(fAnaDebug > 0 ) printf("*** End analysis *** \n");
334   
335 }
336
337 //__________________________________________________________
338 void AliAnaCaloTrackCorrMaker::Terminate(TList * outputList)
339 {  
340   //Execute Terminate of analysis
341   //Do some final plots.
342   
343   if (!outputList) {
344           Error("Terminate", "No output list");
345           return;
346   }
347         
348   for(Int_t iana = 0; iana <  fAnalysisContainer->GetEntries(); iana++){
349     
350     AliAnaCaloTrackCorrBaseClass * ana =  ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
351     if(ana->MakePlotsOn())ana->Terminate(outputList);
352     
353   }//Loop on analysis defined
354 }
355