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