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