]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/CaloTrackCorrBase/AliAnaCaloTrackCorrMaker.cxx
Fix when not being able to open the input file
[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 "TH1F.h"
31 //#include <TObjectTable.h>
32
33 //---- AliRoot system ----
34 #include "AliAnalysisManager.h"
35 #include "AliInputEventHandler.h"
36 #include "AliESDEvent.h"
37 #include "AliAODEvent.h"
38 #include "AliAnaCaloTrackCorrBaseClass.h"
39 #include "AliAnaCaloTrackCorrMaker.h"
40
41 ClassImp(AliAnaCaloTrackCorrMaker)
42
43
44 //__________________________________________________
45 AliAnaCaloTrackCorrMaker::AliAnaCaloTrackCorrMaker() : 
46 TObject(),
47 fReader(0),                   fCaloUtils(0),
48 fOutputContainer(new TList ), fAnalysisContainer(new TList ),
49 fMakeHisto(kFALSE),           fMakeAOD(kFALSE), 
50 fAnaDebug(0),                 fCuts(new TList), 
51 fScaleFactor(-1),
52 fhNEvents(0),                 fhNPileUpEvents(0),
53 fhZVertex(0),                 
54 fhPileUpClusterMult(0),       fhPileUpClusterMultAndSPDPileUp(0),
55 fhTrackMult(0),
56 fhCentrality(0),              fhEventPlaneAngle(0),
57 fhNMergedFiles(0),            fhScaleFactor(0),
58 fhEMCalBCEvent(0),            fhEMCalBCEventCut(0),
59 fhTrackBCEvent(0),            fhTrackBCEventCut(0),
60 fhPrimaryVertexBC(0),         fhTimeStampFraction(0),
61 fhNPileUpVertSPD(0),          fhNPileUpVertTracks(0)
62 {
63   //Default Ctor
64   if(fAnaDebug > 1 ) printf("*** Analysis Maker Constructor *** \n");
65   
66   //Initialize parameters, pointers and histograms
67   InitParameters();
68 }
69
70 //________________________________________________________________________________________
71 AliAnaCaloTrackCorrMaker::AliAnaCaloTrackCorrMaker(const AliAnaCaloTrackCorrMaker & maker) :   
72 TObject(),
73 fReader(),   //(new AliCaloTrackReader(*maker.fReader)),
74 fCaloUtils(),//(new AliCalorimeterUtils(*maker.fCaloUtils)),
75 fOutputContainer(new TList()), fAnalysisContainer(new TList()), 
76 fMakeHisto(maker.fMakeHisto),  fMakeAOD(maker.fMakeAOD),
77 fAnaDebug(maker.fAnaDebug),    fCuts(new TList()),
78 fScaleFactor(maker.fScaleFactor),
79 fhNEvents(maker.fhNEvents), 
80 fhNPileUpEvents(maker.fhNPileUpEvents),
81 fhZVertex(maker.fhZVertex),    
82 fhPileUpClusterMult(maker.fhPileUpClusterMult),
83 fhPileUpClusterMultAndSPDPileUp(maker.fhPileUpClusterMultAndSPDPileUp),
84 fhTrackMult(maker.fhTrackMult),
85 fhCentrality(maker.fhCentrality),
86 fhEventPlaneAngle(maker.fhEventPlaneAngle),
87 fhNMergedFiles(maker.fhNMergedFiles),          
88 fhScaleFactor(maker.fhScaleFactor),
89 fhEMCalBCEvent(maker.fhEMCalBCEvent),
90 fhEMCalBCEventCut(maker.fhEMCalBCEventCut),
91 fhTrackBCEvent(maker.fhTrackBCEvent),
92 fhTrackBCEventCut(maker.fhTrackBCEventCut),
93 fhPrimaryVertexBC(maker.fhPrimaryVertexBC),
94 fhTimeStampFraction(maker.fhTimeStampFraction),
95 fhNPileUpVertSPD(maker.fhNPileUpVertSPD),
96 fhNPileUpVertTracks(maker.fhNPileUpVertTracks)
97 {
98   // cpy ctor
99 }
100
101 //___________________________________________________
102 AliAnaCaloTrackCorrMaker::~AliAnaCaloTrackCorrMaker() 
103 {
104   // Remove all owned pointers.
105   
106   //  Do not delete it here, already done somewhere else, need to understand where.     
107   //  if (fOutputContainer) {
108   //    fOutputContainer->Clear();
109   //    delete fOutputContainer ;
110   //  }   
111   
112   if (fAnalysisContainer)
113   {
114     fAnalysisContainer->Delete();
115     delete fAnalysisContainer ;
116   }   
117   
118   if (fReader)    delete fReader ;
119   if (fCaloUtils) delete fCaloUtils ;
120   
121   if(fCuts)
122   {
123           fCuts->Delete();
124           delete fCuts;
125   }
126         
127 }
128
129 //__________________________________________________________________
130 void    AliAnaCaloTrackCorrMaker::AddAnalysis(TObject* ana, Int_t n) 
131 {
132   // Add analysis depending on AliAnaCaloTrackCorrBaseClass to list
133   
134   if ( fAnalysisContainer)
135   { 
136     fAnalysisContainer->AddAt(ana,n); 
137   }
138   else
139   { 
140     printf("AliAnaCaloTrackCorrMaker::AddAnalysis() - AnalysisContainer not initialized\n");
141     abort();
142   }
143 }  
144
145 //_________________________________________________________
146 TList * AliAnaCaloTrackCorrMaker::FillAndGetAODBranchList()
147
148         
149         // Get any new output AOD branches from analysis and put them in a list
150         // The list is filled in the maker, and new branch passed to the analysis frame
151         // AliAnalysisTaskCaloTrackCorrelation
152   
153         TList *aodBranchList = fReader->GetAODBranchList() ;
154   
155         for(Int_t iana = 0; iana <  fAnalysisContainer->GetEntries(); iana++)
156   {
157                 AliAnaCaloTrackCorrBaseClass * ana =  ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
158                 if(ana->NewOutputAOD()) aodBranchList->Add(ana->GetCreateOutputAODBranch());
159         }
160         
161         return aodBranchList ;
162         
163 }
164
165 //_________________________________________________________
166 void AliAnaCaloTrackCorrMaker::FillControlHistograms()
167 {
168   // Event control histograms
169   
170   AliVEvent* event =  fReader->GetInputEvent();
171   AliESDEvent* esdevent = dynamic_cast<AliESDEvent*> (event);
172   AliAODEvent* aodevent = dynamic_cast<AliAODEvent*> (event);
173   
174   fhNEvents        ->Fill(0); // Number of events analyzed
175   
176   if( fReader->IsPileUpFromSPD())
177     fhNPileUpEvents->Fill(0.5);
178   //if( event->IsPileupFromSPDInMultBins())
179   //  fhNPileUpEvents->Fill(1.5);
180   if( fReader->IsPileUpFromEMCal())
181     fhNPileUpEvents->Fill(2.5);
182   if( fReader->IsPileUpFromSPDOrEMCal() )
183     fhNPileUpEvents->Fill(3.5);
184   if( fReader->IsPileUpFromSPDAndEMCal() )
185     fhNPileUpEvents->Fill(4.5);
186   if( fReader->IsPileUpFromSPDAndNotEMCal() )
187     fhNPileUpEvents->Fill(5.5);
188   if( fReader->IsPileUpFromEMCalAndNotSPD() )
189     fhNPileUpEvents->Fill(6.5);
190   if( fReader->IsPileUpFromNotSPDAndNotEMCal() )
191     fhNPileUpEvents->Fill(7.5);
192   
193   if(fReader->IsPileUpFromSPD())
194     fhPileUpClusterMultAndSPDPileUp ->Fill(fReader->GetNPileUpClusters());
195     
196   fhPileUpClusterMult ->Fill(fReader->GetNPileUpClusters  ());
197   fhTrackMult         ->Fill(fReader->GetTrackMultiplicity());
198   fhCentrality        ->Fill(fReader->GetEventCentrality  ());
199   fhEventPlaneAngle   ->Fill(fReader->GetEventPlaneAngle  ());
200   
201   for(Int_t i = 0; i < 19; i++)
202   {
203     if(fReader->GetTrackEventBC(i))   fhTrackBCEvent   ->Fill(i);
204     if(fReader->GetTrackEventBCcut(i))fhTrackBCEventCut->Fill(i);
205     if(fReader->GetEMCalEventBC(i))   fhEMCalBCEvent   ->Fill(i);
206     if(fReader->GetEMCalEventBCcut(i))fhEMCalBCEventCut->Fill(i);
207   }
208   
209   Double_t v[3];
210   event->GetPrimaryVertex()->GetXYZ(v) ;
211   fhZVertex->Fill(v[2]);
212   
213   Int_t bc = fReader->GetVertexBC();
214   if(bc!=AliVTrack::kTOFBCNA)fhPrimaryVertexBC->Fill(bc);
215   
216   
217   // N pile up vertices
218   Int_t nVerticesSPD    = -1;
219   Int_t nVerticesTracks = -1;
220   
221   if      (esdevent)
222   {
223     nVerticesSPD    = esdevent->GetNumberOfPileupVerticesSPD();
224     nVerticesTracks = esdevent->GetNumberOfPileupVerticesTracks();
225     
226   }//ESD
227   else if (aodevent)
228   {
229     nVerticesSPD    = aodevent->GetNumberOfPileupVerticesSPD();
230     nVerticesTracks = aodevent->GetNumberOfPileupVerticesTracks();
231   }//AOD
232   
233   fhNPileUpVertSPD   ->Fill(nVerticesSPD);
234   fhNPileUpVertTracks->Fill(nVerticesTracks);
235
236   // Time stamp
237   if(fReader->IsSelectEventTimeStampOn() && esdevent)
238   {
239     Int_t timeStamp = esdevent->GetTimeStamp();
240     Float_t timeStampFrac = 1.*(timeStamp-fReader->GetRunTimeStampMin()) /
241                                (fReader->GetRunTimeStampMax()-fReader->GetRunTimeStampMin());
242     
243     //printf("stamp %d, min %d, max %d, frac %f\n", timeStamp, fReader->GetRunTimeStampMin(), fReader->GetRunTimeStampMax(), timeStampFrac);
244
245     fhTimeStampFraction->Fill(timeStampFrac);
246   }
247   
248   
249 }
250
251 //_______________________________________________________
252 TList * AliAnaCaloTrackCorrMaker::GetListOfAnalysisCuts()
253
254   
255         // Get the list of the cuts used for the analysis
256         // The list is filled in the maker, called by the task in LocalInit() and posted there
257   
258         for(Int_t iana = 0; iana <  fAnalysisContainer->GetEntries(); iana++)
259   {
260                 AliAnaCaloTrackCorrBaseClass * ana =  ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
261                 TObjString * objstring = ana->GetAnalysisCuts();
262     
263                 if(objstring)fCuts->Add(objstring);
264         }
265   
266         return fCuts ;
267   
268 }
269
270 //___________________________________________________
271 TList *AliAnaCaloTrackCorrMaker::GetOutputContainer()
272 {
273   // Fill the output list of histograms during the CreateOutputObjects stage.
274   
275   //Initialize calorimeters  geometry pointers
276   //GetCaloUtils()->InitPHOSGeometry();
277   //GetCaloUtils()->InitEMCALGeometry();
278   
279   //General event histograms
280   
281   fhNEvents      = new TH1F("hNEvents",   "Number of analyzed events"     , 1 , 0 , 1  ) ;
282   fhNEvents->SetYTitle("# events");
283   fOutputContainer->Add(fhNEvents);
284   
285   fhNPileUpEvents      = new TH1F("hNPileUpEvents",   "Number of events considered as pile-up", 8 , 0 , 8 ) ;
286   fhNPileUpEvents->SetYTitle("# events");
287   fhNPileUpEvents->GetXaxis()->SetBinLabel(1 ,"SPD");
288   fhNPileUpEvents->GetXaxis()->SetBinLabel(2 ,"Multi SPD");
289   fhNPileUpEvents->GetXaxis()->SetBinLabel(3 ,"EMCal");
290   fhNPileUpEvents->GetXaxis()->SetBinLabel(4 ,"EMCal || SPD");
291   fhNPileUpEvents->GetXaxis()->SetBinLabel(5 ,"EMCal && SPD");
292   fhNPileUpEvents->GetXaxis()->SetBinLabel(6 ,"!EMCal && SPD");
293   fhNPileUpEvents->GetXaxis()->SetBinLabel(7 ,"EMCal && !SPD");
294   fhNPileUpEvents->GetXaxis()->SetBinLabel(8 ,"!EMCal && !SPD");
295   fOutputContainer->Add(fhNPileUpEvents);
296   
297   fhTrackBCEvent      = new TH1F("hTrackBCEvent",   "Number of events with at least 1 track in a bunch crossing ", 19 , 0 , 19 ) ;
298   fhTrackBCEvent->SetYTitle("# events");
299   fhTrackBCEvent->SetXTitle("Bunch crossing");
300   for(Int_t i = 1; i < 20; i++)
301     fhTrackBCEvent->GetXaxis()->SetBinLabel(i ,Form("%d",i-10));
302   fOutputContainer->Add(fhTrackBCEvent);
303   
304   fhTrackBCEventCut      = new TH1F("hTrackBCEventCut",   "Number of events with at least 1 track in a bunch crossing ", 19 , 0 , 19 ) ;
305   fhTrackBCEventCut->SetYTitle("# events");
306   fhTrackBCEventCut->SetXTitle("Bunch crossing");
307   for(Int_t i = 1; i < 20; i++)
308     fhTrackBCEventCut->GetXaxis()->SetBinLabel(i ,Form("%d",i-10));
309   fOutputContainer->Add(fhTrackBCEventCut);
310
311   fhPrimaryVertexBC      = new TH1F("hPrimaryVertexBC", "Number of primary vertex per bunch crossing ", 41 , -20 , 20  ) ;
312   fhPrimaryVertexBC->SetYTitle("# events");
313   fhPrimaryVertexBC->SetXTitle("Bunch crossing");
314   fOutputContainer->Add(fhPrimaryVertexBC);
315
316   fhEMCalBCEvent      = new TH1F("hEMCalBCEvent",   "Number of events with at least 1 cluster in a bunch crossing ", 19 , 0 , 19 ) ;
317   fhEMCalBCEvent->SetYTitle("# events");
318   fhEMCalBCEvent->SetXTitle("Bunch crossing");
319   for(Int_t i = 1; i < 20; i++)
320     fhEMCalBCEvent->GetXaxis()->SetBinLabel(i ,Form("%d",i-10));
321   fOutputContainer->Add(fhEMCalBCEvent);
322   
323   fhEMCalBCEventCut      = new TH1F("hEMCalBCEventCut",   "Number of events with at least 1 cluster in a bunch crossing", 19 , 0 , 19 ) ;
324   fhEMCalBCEventCut->SetYTitle("# events");
325   fhEMCalBCEventCut->SetXTitle("Bunch crossing");
326   for(Int_t i = 1; i < 20; i++)
327     fhEMCalBCEventCut->GetXaxis()->SetBinLabel(i ,Form("%d",i-10));
328   fOutputContainer->Add(fhEMCalBCEventCut);
329   
330   fhZVertex      = new TH1F("hZVertex", " Z vertex distribution"   , 200 , -50 , 50  ) ;
331   fhZVertex->SetXTitle("v_{z} (cm)");
332   fOutputContainer->Add(fhZVertex);
333   
334   fhTrackMult    = new TH1F("hTrackMult", "Number of tracks per events"   , 2000 , 0 , 2000  ) ;
335   fhTrackMult->SetXTitle("# tracks");
336   fOutputContainer->Add(fhTrackMult);
337   
338   fhPileUpClusterMult    = new TH1F("hPileUpClusterMult", "Number of clusters per event with large time (|t| > 20 ns)" , 100 , 0 , 100  ) ;
339   fhPileUpClusterMult->SetXTitle("# clusters");
340   fOutputContainer->Add(fhPileUpClusterMult);
341   
342   fhPileUpClusterMultAndSPDPileUp = new TH1F("hPileUpClusterMultAndSPDPileUp", "Number of clusters per event with large time (|t| > 20 ns, events tagged as pile-up by SPD)" , 100 , 0 , 100 ) ;
343   fhPileUpClusterMultAndSPDPileUp->SetXTitle("# clusters");
344   fOutputContainer->Add(fhPileUpClusterMultAndSPDPileUp);
345   
346   fhNPileUpVertSPD  = new TH1F ("hNPileUpVertSPD","N pile-up SPD vertex", 50,0,50);
347   fhNPileUpVertSPD->SetYTitle("# vertex ");
348   fOutputContainer->Add(fhNPileUpVertSPD);
349   
350   fhNPileUpVertTracks  = new TH1F ("hNPileUpVertTracks","N pile-up Tracks vertex", 50,0,50);
351   fhNPileUpVertTracks->SetYTitle("# vertex ");
352   fOutputContainer->Add(fhNPileUpVertTracks);
353   
354   fhCentrality   = new TH1F("hCentrality","Number of events in centrality bin",100,0.,100) ;
355   fhCentrality->SetXTitle("Centrality bin");
356   fOutputContainer->Add(fhCentrality) ;  
357   
358   fhEventPlaneAngle=new TH1F("hEventPlaneAngle","Number of events in event plane",100,0.,TMath::Pi()) ;
359   fhEventPlaneAngle->SetXTitle("EP angle (rad)");
360   fOutputContainer->Add(fhEventPlaneAngle) ;
361   
362   if(fReader->IsSelectEventTimeStampOn())
363   {
364     fhTimeStampFraction = new TH1F("hTimeStampFraction","Fraction of events within a given time stamp range",150, -1, 2) ;
365     fhTimeStampFraction->SetXTitle("fraction");
366     fOutputContainer->Add(fhTimeStampFraction) ;
367   }
368   
369   if(fScaleFactor > 0)
370   {
371     fhNMergedFiles = new TH1F("hNMergedFiles",   "Number of merged output files"     , 1 , 0 , 1  ) ;
372     fhNMergedFiles->SetYTitle("# files");
373     fhNMergedFiles->Fill(1); // Fill here with one entry, while merging it will count the rest
374     fOutputContainer->Add(fhNMergedFiles);
375     
376     fhScaleFactor = new TH1F("hScaleFactor",   "Number of merged output files"     , 1 , 0 , 1  ) ;
377     fhScaleFactor->SetYTitle("scale factor");
378     fhScaleFactor->SetBinContent(1,fScaleFactor); // Fill here 
379     fOutputContainer->Add(fhScaleFactor);    
380   }
381   
382   if(!fAnalysisContainer || fAnalysisContainer->GetEntries()==0)
383   {
384     printf("AliAnaCaloTrackCorrMaker::GetOutputContainer() - Analysis job list not initialized!!!\n");
385     return fOutputContainer;
386   }
387   
388   const Int_t buffersize = 255;
389   char newname[buffersize];
390   for(Int_t iana = 0; iana <  fAnalysisContainer->GetEntries(); iana++)
391   {
392     
393     AliAnaCaloTrackCorrBaseClass * ana =  ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
394     
395     if(fMakeHisto) // Analysis with histograms as output on
396     {
397       
398       //Fill container with appropriate histograms                      
399       TList * templist =  ana ->GetCreateOutputObjects(); 
400       templist->SetOwner(kFALSE); //Owner is fOutputContainer.
401       
402       for(Int_t i = 0; i < templist->GetEntries(); i++)
403       {
404         
405         //Add only  to the histogram name the name of the task
406         if(   strcmp((templist->At(i))->ClassName(),"TObjString")   ) 
407         {
408           snprintf(newname,buffersize, "%s%s", (ana->GetAddedHistogramsStringToName()).Data(), (templist->At(i))->GetName());
409           //printf("name %s, new name %s\n",(templist->At(i))->GetName(),newname);
410           ((TH1*) templist->At(i))->SetName(newname);
411         }
412         
413         //Add histogram to general container
414         fOutputContainer->Add(templist->At(i)) ;
415         
416       }
417       
418       delete templist;
419       
420     }// Analysis with histograms as output on
421     
422   }//Loop on analysis defined
423   
424   return fOutputContainer;
425   
426 }
427
428 //___________________________________
429 void AliAnaCaloTrackCorrMaker::Init()
430 {  
431   //Init container histograms and other common variables
432   // Fill the output list of histograms during the CreateOutputObjects stage.
433   
434   //Initialize reader
435   GetReader()->Init();
436   GetReader()->SetCaloUtils(GetCaloUtils()); // pass the calo utils pointer to the reader
437         
438   
439   if(!fAnalysisContainer || fAnalysisContainer->GetEntries()==0)
440   {
441     printf("AliAnaCaloTrackCorrMaker::GetOutputInit() - Analysis job list not initialized!!!\n");
442     return;
443   }
444   
445   for(Int_t iana = 0; iana <  fAnalysisContainer->GetEntries(); iana++)
446   {
447     
448     AliAnaCaloTrackCorrBaseClass * ana =  ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
449     
450     ana->SetReader(fReader);       // Set Reader for each analysis
451     ana->SetCaloUtils(fCaloUtils); // Set CaloUtils for each analysis
452     
453     ana->Init();
454     
455   }//Loop on analysis defined
456   
457 }
458
459 //_____________________________________________
460 void AliAnaCaloTrackCorrMaker::InitParameters()
461 {       
462   //Init data members
463   
464   fMakeHisto  = kTRUE;
465   fMakeAOD    = kTRUE; 
466   fAnaDebug   = 0; // No debugging info displayed by default
467         
468 }
469
470 //______________________________________________________________
471 void AliAnaCaloTrackCorrMaker::Print(const Option_t * opt) const
472 {       
473   //Print some relevant parameters set for the analysis
474         
475   if(! opt)
476     return;
477   
478   printf("***** Print: %s %s ******\n", GetName(), GetTitle() ) ;
479   printf("Debug level                =     %d\n", fAnaDebug   ) ;
480   printf("Produce Histo              =     %d\n", fMakeHisto  ) ;
481   printf("Produce AOD                =     %d\n", fMakeAOD    ) ;
482   printf("Number of analysis tasks   =     %d\n", fAnalysisContainer->GetEntries()) ;
483   
484   if(!strcmp("all",opt))
485   {
486           printf("Print analysis Tasks settings :\n") ;
487           for(Int_t iana = 0; iana<fAnalysisContainer->GetEntries(); iana++)
488     {
489                   ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana))->Print("");
490           }
491     
492           printf("Print analysis Reader settings :\n") ;
493           fReader->Print("");
494           printf("Print analysis Calorimeter Utils settings :\n") ;
495           fCaloUtils->Print("");
496     
497   }
498   
499
500
501 //_______________________________________________________________________
502 void AliAnaCaloTrackCorrMaker::ProcessEvent(const Int_t iEntry, 
503                                             const char * currentFileName)
504 {
505   //Process analysis for this event
506   
507   if(fMakeHisto && !fOutputContainer)
508   {
509     printf("AliAnaCaloTrackCorrMaker::ProcessEvent() - Histograms not initialized\n");
510     abort();
511   }
512         
513   if(fAnaDebug >= 0 )
514   {
515                 printf("***  AliAnaCaloTrackCorrMaker::ProcessEvent() Event %d   ***  \n",iEntry);
516           if(fAnaDebug > 1 ) 
517     {
518                   printf("AliAnaCaloTrackCorrMaker::ProcessEvent() - Current File Name : %s\n", currentFileName);
519                   //printf("fAODBranchList %p, entries %d\n",fAODBranchList,fAODBranchList->GetEntries());
520           }
521   }
522   
523   //Each event needs an empty branch
524   TList * aodList = fReader->GetAODBranchList();
525   Int_t nAODBranches = aodList->GetEntries();
526   for(Int_t iaod = 0; iaod < nAODBranches; iaod++)
527   {
528           TClonesArray *tca = dynamic_cast<TClonesArray*> (aodList->At(iaod));
529           if(tca) tca->Clear("C");
530   }
531   
532   //Set geometry matrices before filling arrays, in case recalibration/position calculation etc is needed
533   fCaloUtils->AccessGeometry(fReader->GetInputEvent()); 
534   
535   //Set the AODB calibration, bad channels etc. parameters at least once
536   fCaloUtils->AccessOADB(fReader->GetInputEvent());     
537
538   
539   //Tell the reader to fill the data in the 3 detector lists
540   Bool_t ok = fReader->FillInputEvent(iEntry, currentFileName);
541   if(!ok)
542   {
543           if(fAnaDebug >= 1 )printf("*** Skip event *** %d \n",iEntry);
544     fReader->ResetLists();
545           return ;
546   }
547         
548   //Magic line to write events to file
549   if(fReader->WriteDeltaAODToFile())AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler()->SetFillAOD(kTRUE);
550   
551   //printf(">>>>>>>>>> BEFORE >>>>>>>>>>>\n");
552   //gObjectTable->Print();
553   
554   //Access pointers, and trigger mask check needed in mixing case
555   AliAnalysisManager   *manager      = AliAnalysisManager::GetAnalysisManager();
556   AliInputEventHandler *inputHandler = dynamic_cast<AliInputEventHandler*>(manager->GetInputEventHandler());
557   
558   UInt_t isMBTrigger = kFALSE;
559   UInt_t isTrigger   = kFALSE;
560   if(inputHandler)
561   {
562     isMBTrigger = inputHandler->IsEventSelected() & fReader->GetMixEventTriggerMask();
563     isTrigger   = inputHandler->IsEventSelected() & fReader->GetEventTriggerMask();
564   }
565   
566   //Loop on analysis algorithms
567   
568   if(fAnaDebug > 0 ) printf("*** Begin analysis *** \n");
569   
570   Int_t nana = fAnalysisContainer->GetEntries() ;
571   for(Int_t iana = 0; iana <  nana; iana++)
572   {
573     AliAnaCaloTrackCorrBaseClass * ana =  ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ; 
574     
575     ana->ConnectInputOutputAODBranches(); //Sets branches for each analysis
576     
577     //Fill pool for mixed event for the analysis that need it
578     if(!fReader->IsEventTriggerAtSEOn() && isMBTrigger)
579     {
580       ana->FillEventMixPool();
581       continue; // pool filled do not try to fill AODs or histograms
582     }
583     
584     //Make analysis, create aods in aod branch and in some cases fill histograms
585     if(fMakeAOD  )  ana->MakeAnalysisFillAOD()  ;
586     
587     //Make further analysis with aod branch and fill histograms
588     if(fMakeHisto)  ana->MakeAnalysisFillHistograms()  ;
589     
590   }
591         
592   fReader->ResetLists();
593
594   // In case of mixing analysis, non triggered events are used,
595   // do not fill control histograms for a non requested triggered event
596   if(!fReader->IsEventTriggerAtSEOn() && !isTrigger)
597   {    
598     if(fAnaDebug > 0 ) printf("AliAnaCaloTrackMaker::ProcessEvent() - *** End analysis, MB for mixing *** \n");
599     return;
600   }
601   
602   FillControlHistograms();
603   
604   //printf(">>>>>>>>>> AFTER >>>>>>>>>>>\n");
605   //gObjectTable->Print();
606         
607   if(fAnaDebug > 0 ) printf("AliAnaCaloTrackMaker::ProcessEvent() - *** End analysis *** \n");
608   
609 }
610
611 //__________________________________________________________
612 void AliAnaCaloTrackCorrMaker::Terminate(TList * outputList)
613 {  
614   //Execute Terminate of analysis
615   //Do some final plots.
616   
617   if (!outputList) 
618   {
619           Error("Terminate", "No output list");
620           return;
621   }
622           
623   for(Int_t iana = 0; iana <  fAnalysisContainer->GetEntries(); iana++)
624   {
625     
626     AliAnaCaloTrackCorrBaseClass * ana =  ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
627     if(ana->MakePlotsOn())ana->Terminate(outputList);
628     
629   }//Loop on analysis defined
630   
631 }
632