]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/PartCorrBase/AliCaloTrackReader.cxx
0d1caae08f8f2b55853aa2fc7b09e4baba1e8ed1
[u/mrichter/AliRoot.git] / PWG4 / PartCorrBase / AliCaloTrackReader.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 // Base class for reading data: MonteCarlo, ESD or AOD, of PHOS EMCAL and 
19 // Central Barrel Tracking detectors (CTS).
20 // Not all MC particles/tracks/clusters are kept, some kinematical/fiducial restrictions are done.
21 // Mother class of : AliCaloTrackESDReader: Fills ESD data in 3 TObjArrays (PHOS, EMCAL, CTS)
22 //                 : AliCaloTrackMCReader: Fills Kinematics data in 3 TObjArrays (PHOS, EMCAL, CTS)
23 //                 : AliCaloTrackAODReader: Fills AOD data in 3 TObjArrays (PHOS, EMCAL, CTS) 
24 //                
25 //-- Author: Gustavo Conesa (LNF-INFN) 
26 //////////////////////////////////////////////////////////////////////////////
27
28
29 // --- ROOT system ---
30 #include "TFile.h"
31
32 //---- ANALYSIS system ----
33 #include "AliCaloTrackReader.h"
34 #include "AliMCEvent.h"
35 #include "AliAODMCHeader.h"
36 #include "AliGenPythiaEventHeader.h"
37 #include "AliAODEvent.h"
38
39 ClassImp(AliCaloTrackReader)
40   
41   
42 //____________________________________________________________________________
43   AliCaloTrackReader::AliCaloTrackReader() : 
44     TObject(), fEventNumber(-1), fCurrentFileName(""),fDataType(0), fDebug(0), 
45     fFiducialCut(0x0), fCheckFidCut(kFALSE), fComparePtHardAndJetPt(kFALSE), fPtHardAndJetPtFactor(7),
46     fCTSPtMin(0), fEMCALPtMin(0),fPHOSPtMin(0),
47     fAODCTS(new TObjArray()), fAODEMCAL(new TObjArray()), fAODPHOS(new TObjArray()),
48     fEMCALCells(0x0), fPHOSCells(0x0),
49     fInputEvent(0x0), fOutputEvent(0x0),fMC(0x0),
50     fFillCTS(0),fFillEMCAL(0),fFillPHOS(0),
51     fFillEMCALCells(0),fFillPHOSCells(0), 
52     fSecondInputAODTree(0x0), fSecondInputAODEvent(0x0),
53     fSecondInputFileName(""),fSecondInputFirstEvent(0), 
54     fAODCTSNormalInputEntries(0), fAODEMCALNormalInputEntries(0), 
55     fAODPHOSNormalInputEntries(0), fTrackStatus(0), 
56     fReadStack(kFALSE), fReadAODMCParticles(kFALSE), 
57     fCleanOutputStdAOD(kFALSE), fDeltaAODFileName("deltaAODPartCorr.root"),fFiredTriggerClassName(""),
58     fAnaLED(kFALSE),fTaskName(""),fCaloUtils(0x0)
59 {
60   //Ctor
61   
62   //Initialize parameters
63   InitParameters();
64 }
65 /*
66 //____________________________________________________________________________
67 AliCaloTrackReader::AliCaloTrackReader(const AliCaloTrackReader & reader) :   
68   TObject(reader), fEventNumber(reader.fEventNumber), fCurrentFileName(reader.fCurrentFileName), 
69   fDataType(reader.fDataType), fDebug(reader.fDebug),
70   fFiducialCut(reader.fFiducialCut),
71   fComparePtHardAndJetPt(reader.fComparePtHardAndJetPt),
72   fPtHardAndJetPtFactor(reader.fPtHardAndJetPtFactor),
73   fCTSPtMin(reader.fCTSPtMin), fEMCALPtMin(reader.fEMCALPtMin),fPHOSPtMin(reader.fPHOSPtMin), 
74   fAODCTS(new TObjArray(*reader.fAODCTS)),  
75   fAODEMCAL(new TObjArray(*reader.fAODEMCAL)),
76   fAODPHOS(new TObjArray(*reader.fAODPHOS)),
77   fEMCALCells(new TNamed(*reader.fEMCALCells)),
78   fPHOSCells(new TNamed(*reader.fPHOSCells)),
79   fInputEvent(reader.fInputEvent), fOutputEvent(reader.fOutputEvent), fMC(reader.fMC),
80   fFillCTS(reader.fFillCTS),fFillEMCAL(reader.fFillEMCAL),fFillPHOS(reader.fFillPHOS),
81   fFillEMCALCells(reader.fFillEMCALCells),fFillPHOSCells(reader.fFillPHOSCells),
82   fSecondInputAODTree(reader.fSecondInputAODTree), 
83   fSecondInputAODEvent(reader.fSecondInputAODEvent),
84   fSecondInputFileName(reader.fSecondInputFileName), 
85   fSecondInputFirstEvent(reader.fSecondInputFirstEvent),
86   fAODCTSNormalInputEntries(reader.fAODCTSNormalInputEntries), 
87   fAODEMCALNormalInputEntries(reader.fAODEMCALNormalInputEntries), 
88   fAODPHOSNormalInputEntries(reader.fAODPHOSNormalInputEntries),
89   fTrackStatus(reader.fTrackStatus),
90   fReadStack(reader.fReadStack), fReadAODMCParticles(reader.fReadAODMCParticles),
91   fCleanOutputStdAOD(reader.fCleanOutputStdAOD), fDeltaAODFileName(reader.fDeltaAODFileName),
92   fFiredTriggerClassName(reader.fFiredTriggerClassName),
93   fAnaLED(reader.fAnaLED), 
94   fTaskName(reader.fTaskName),
95   fCaloUtils(new AliCalorimeterUtils(*reader.fCaloUtils))
96 {
97   // cpy ctor  
98 }
99 */
100 //_________________________________________________________________________
101 //AliCaloTrackReader & AliCaloTrackReader::operator = (const AliCaloTrackReader & source)
102 //{
103 //  // assignment operator
104 //  
105 //  if(&source == this) return *this;
106 //  
107 //  fDataType    = source.fDataType ;
108 //  fDebug       = source.fDebug ;
109 //  fEventNumber = source.fEventNumber ;
110 //  fCurrentFileName = source.fCurrentFileName ;
111 //  fFiducialCut = source.fFiducialCut;
112 //      
113 //  fComparePtHardAndJetPt = source.fComparePtHardAndJetPt;
114 //  fPtHardAndJetPtFactor  = source.fPtHardAndJetPtFactor;
115 //      
116 //  fCTSPtMin    = source.fCTSPtMin ;
117 //  fEMCALPtMin  = source.fEMCALPtMin ;
118 //  fPHOSPtMin   = source.fPHOSPtMin ; 
119 //  
120 //  fAODCTS     = new TObjArray(*source.fAODCTS) ;
121 //  fAODEMCAL   = new TObjArray(*source.fAODEMCAL) ;
122 //  fAODPHOS    = new TObjArray(*source.fAODPHOS) ;
123 //  fEMCALCells = new TNamed(*source.fEMCALCells) ;
124 //  fPHOSCells  = new TNamed(*source.fPHOSCells) ;
125 //
126 //  fInputEvent  = source.fInputEvent;
127 //  fOutputEvent = source.fOutputEvent;
128 //  fMC          = source.fMC;
129 //  
130 //  fFillCTS        = source.fFillCTS;
131 //  fFillEMCAL      = source.fFillEMCAL;
132 //  fFillPHOS       = source.fFillPHOS;
133 //  fFillEMCALCells = source.fFillEMCALCells;
134 //  fFillPHOSCells  = source.fFillPHOSCells;
135 //
136 //  fSecondInputAODTree    = source.fSecondInputAODTree;
137 //  fSecondInputAODEvent   = source.fSecondInputAODEvent;
138 //  fSecondInputFileName   = source.fSecondInputFileName;
139 //  fSecondInputFirstEvent = source.fSecondInputFirstEvent;
140 //
141 //  fAODCTSNormalInputEntries   = source.fAODCTSNormalInputEntries; 
142 //  fAODEMCALNormalInputEntries = source.fAODEMCALNormalInputEntries; 
143 //  fAODPHOSNormalInputEntries  = source.fAODPHOSNormalInputEntries;
144 //      
145 //  fTrackStatus        = source.fTrackStatus;
146 //  fReadStack          = source.fReadStack;
147 //  fReadAODMCParticles = source.fReadAODMCParticles;   
148 //      
149 //  fCleanOutputStdAOD  = source.fCleanOutputStdAOD;
150 //  fDeltaAODFileName   = source.fDeltaAODFileName;
151 //      
152 //  fFiredTriggerClassName = source.fFiredTriggerClassName  ;
153 //      
154 //  return *this;
155 //  
156 //}
157
158 //_________________________________
159 AliCaloTrackReader::~AliCaloTrackReader() {
160   //Dtor
161   
162   if(fFiducialCut) delete fFiducialCut ;
163         
164   if(fAODCTS){
165     fAODCTS->Clear() ; 
166     delete fAODCTS ;
167   }
168   
169   if(fAODEMCAL){
170     fAODEMCAL->Clear() ; 
171     delete fAODEMCAL ;
172   }
173   
174   if(fAODPHOS){
175     fAODPHOS->Clear() ; 
176     delete fAODPHOS ;
177   }
178   
179   if(fEMCALCells){
180     fEMCALCells->Clear() ; 
181     delete fEMCALCells ;
182   }
183   
184   if(fPHOSCells){
185     fPHOSCells->Clear() ; 
186     delete fPHOSCells ;
187   }
188         
189 //  Pointers not owned, done by the analysis frame
190 //  if(fInputEvent)  delete fInputEvent ;
191 //  if(fOutputEvent) delete fOutputEvent ;
192 //  if(fMC)          delete fMC ;  
193         
194 //  if(fSecondInputAODTree){
195 //    fSecondInputAODTree->Clear();
196 //    delete fSecondInputAODTree;
197 //  }
198 //      
199 //  if(fSecondInputAODEvent) delete fSecondInputAODEvent ;
200         
201   //  Pointer not owned, deleted by maker
202   //if (fCaloUtils) delete fCaloUtils ;
203
204 }
205
206 //_________________________________________________________________________
207 Bool_t AliCaloTrackReader::ComparePtHardAndJetPt(){
208         // Check the event, if the requested ptHard is much larger than the jet pT, then there is a problem.
209         // Only for PYTHIA.
210         if(!fReadStack) return kTRUE; //Information not filtered to AOD
211         
212         if(!strcmp(GetGenEventHeader()->ClassName(), "AliGenPythiaEventHeader")){
213                 TParticle * jet =  new TParticle;
214                 AliGenPythiaEventHeader* pygeh= (AliGenPythiaEventHeader*) GetGenEventHeader();
215                 Int_t nTriggerJets =  pygeh->NTriggerJets();
216                 Float_t ptHard = pygeh->GetPtHard();
217                 
218                 //if(fDebug > 1) printf("AliMCAnalysisUtils::PythiaEventHeader: Njets: %d, pT Hard %f\n",nTriggerJets, ptHard);
219             Float_t tmpjet[]={0,0,0,0};
220                 for(Int_t ijet = 0; ijet< nTriggerJets; ijet++){
221                         pygeh->TriggerJet(ijet, tmpjet);
222                         jet = new TParticle(94, 21, -1, -1, -1, -1, tmpjet[0],tmpjet[1],tmpjet[2],tmpjet[3], 0,0,0,0);
223                         //Compare jet pT and pt Hard
224                         //if(fDebug > 1) printf("AliMCAnalysisUtils:: %d pycell jet pT %f\n",ijet, jet->Pt());
225                         if(jet->Pt() > fPtHardAndJetPtFactor * ptHard) {
226                                 printf("AliMCAnalysisUtils::PythiaEventHeader: Njets: %d, pT Hard %2.2f, pycell jet pT %2.2f, rejection factor %1.1f\n",
227                                            nTriggerJets, ptHard, jet->Pt(), fPtHardAndJetPtFactor);
228                                 return kFALSE;
229                         }
230                 }
231         }
232         
233         return kTRUE ;
234         
235 }
236
237 //____________________________________________________________________________
238 AliStack* AliCaloTrackReader::GetStack() const {
239   //Return pointer to stack
240   if(fMC)
241     return fMC->Stack();
242   else{
243     if(fDebug > 1) printf("AliCaloTrackReader::GetStack() - Stack is not available\n"); 
244     return 0x0 ;
245   }
246 }
247
248 //____________________________________________________________________________
249 AliHeader* AliCaloTrackReader::GetHeader() const {
250   //Return pointer to header
251   if(fMC)
252     return fMC->Header();
253   else{
254     printf("AliCaloTrackReader::Header is not available\n"); 
255     return 0x0 ;
256   }
257 }
258 //____________________________________________________________________________
259 AliGenEventHeader* AliCaloTrackReader::GetGenEventHeader() const {
260   //Return pointer to Generated event header
261   if(fMC)
262     return fMC->GenEventHeader();
263   else{
264     printf("AliCaloTrackReader::GenEventHeader is not available\n"); 
265     return 0x0 ;
266   }
267 }
268
269 //____________________________________________________________________________
270 TClonesArray* AliCaloTrackReader::GetAODMCParticles(Int_t input) const {
271         //Return list of particles in AOD. Do it for the corresponding input event.
272         if(fDataType == kAOD){
273          //Normal input AOD
274          if(input == 0) return (TClonesArray*)((AliAODEvent*)fInputEvent)->FindListObject("mcparticles");
275           //Second input AOD
276          else if(input == 1 && fSecondInputAODEvent) return (TClonesArray*) fSecondInputAODEvent->FindListObject("mcparticles");        
277          else {
278              printf("AliCaloTrackReader::GetAODMCParticles() - wrong AOD input index? %d, or non existing tree? \n",input); 
279                  return 0x0;
280          }
281         }
282         else {
283                 printf("AliCaloTrackReader::GetAODMCParticles() - Input are not AODs\n"); 
284                 return 0x0;
285         }
286 }
287
288 //____________________________________________________________________________
289 AliAODMCHeader* AliCaloTrackReader::GetAODMCHeader(Int_t input) const {
290         //Return MC header in AOD. Do it for the corresponding input event.
291         if(fDataType == kAOD){
292                 //Normal input AOD
293                 if(input == 0) return (AliAODMCHeader*)((AliAODEvent*)fInputEvent)->FindListObject("mcheader");
294                 //Second input AOD
295                 else if(input == 1) return  (AliAODMCHeader*) fSecondInputAODEvent->FindListObject("mcheader"); 
296                 else {
297                         printf("AliCaloTrackReader::GetAODMCHeader() - wrong AOD input index, %d\n",input);
298                         return 0x0;
299                 }
300         }
301         else {
302                 printf("AliCaloTrackReader::GetAODMCHeader() - Input are not AODs\n");
303                 return 0x0;
304         }
305 }
306
307 //_______________________________________________________________
308 void AliCaloTrackReader::Init()
309 {
310         //Init reader. Method to be called in AliAnaPartCorrMaker
311         
312         //Get the file with second input events if the filename is given
313         //Get the tree and connect the AODEvent. Only with AODs
314
315         if(fReadStack && fReadAODMCParticles){
316                 printf("AliCaloTrackReader::Init() - Cannot access stack and mcparticles at the same time, change them \n");
317                 fReadStack = kFALSE;
318                 fReadAODMCParticles = kFALSE;
319         }
320         
321         if(fSecondInputFileName!=""){
322                 if(fDataType == kAOD){
323                         TFile * input2   = new TFile(fSecondInputFileName,"read");
324                         printf("AliCaloTrackReader::Init() - Second input file opened: %s, size %d \n", input2->GetName(), (Int_t) input2->GetSize());
325                         fSecondInputAODTree = (TTree*) input2->Get("aodTree");
326                         if(fSecondInputAODTree) printf("AliCaloTrackReader::Init() - Second input tree opened: %s, entries %d \n", 
327                                                                                    fSecondInputAODTree->GetName(), (Int_t) fSecondInputAODTree->GetEntries());
328                         else{
329                          printf("AliCaloTrackReader::Init() - Second input tree not available, STOP \n");
330                          abort();
331                         }
332                         fSecondInputAODEvent = new AliAODEvent;
333                         fSecondInputAODEvent->ReadFromTree(fSecondInputAODTree);
334                         if(fSecondInputFirstEvent >= fSecondInputAODTree->GetEntriesFast()){
335                                 printf("AliCaloTrackReader::Init() - Requested first event of second input %d, is larger than number of events %d, STOP\n", 
336                                            fSecondInputFirstEvent, (Int_t) fSecondInputAODTree->GetEntriesFast());
337                                 abort();
338                         }
339                 }
340                 else printf("AliCaloTrackReader::Init() - Second input not added, reader is not AOD\n");
341         }
342 }
343
344 //_______________________________________________________________
345 void AliCaloTrackReader::InitParameters()
346 {
347   //Initialize the parameters of the analysis.
348   fDataType   = kESD ;
349   fCTSPtMin   = 0.2 ;
350   fEMCALPtMin = 0.2 ;
351   fPHOSPtMin  = 0.2 ;
352
353   //Do not filter the detectors input by default.
354   fFillEMCAL      = kFALSE;
355   fFillPHOS       = kFALSE;
356   fFillCTS        = kFALSE;
357   fFillEMCALCells = kFALSE;
358   fFillPHOSCells  = kFALSE;
359
360   fSecondInputFileName   = "" ;
361   fSecondInputFirstEvent = 0 ;
362   fReadStack             = kFALSE; // Check in the constructor of the other readers if it was set or in the configuration file
363   fReadAODMCParticles    = kFALSE; // Check in the constructor of the other readers if it was set or in the configuration file
364   fCleanOutputStdAOD     = kFALSE; // Clean the standard clusters/tracks?
365   fDeltaAODFileName      = "deltaAODPartCorr.root";
366   fFiredTriggerClassName      = "";
367                 
368   fAnaLED = kFALSE;
369 }
370
371 //________________________________________________________________
372 void AliCaloTrackReader::Print(const Option_t * opt) const
373 {
374
375   //Print some relevant parameters set for the analysis
376   if(! opt)
377     return;
378
379   printf("***** Print: %s %s ******\n", GetName(), GetTitle() ) ;
380   printf("Task name      : %s\n", fTaskName.Data()) ;
381   printf("Data type      : %d\n", fDataType) ;
382   printf("CTS Min pT     : %2.1f GeV/c\n", fCTSPtMin) ;
383   printf("EMCAL Min pT   : %2.1f GeV/c\n", fEMCALPtMin) ;
384   printf("PHOS Min pT    : %2.1f GeV/c\n", fPHOSPtMin) ;
385   printf("Use CTS         =     %d\n", fFillCTS) ;
386   printf("Use EMCAL       =     %d\n", fFillEMCAL) ;
387   printf("Use PHOS        =     %d\n", fFillPHOS) ;
388   printf("Use EMCAL Cells =     %d\n", fFillEMCALCells) ;
389   printf("Use PHOS  Cells =     %d\n", fFillPHOSCells) ;
390   printf("Track status    =     %d\n", (Int_t) fTrackStatus) ;
391   if(fComparePtHardAndJetPt)
392           printf("Compare jet pt and pt hard to accept event, factor = %2.2f",fPtHardAndJetPtFactor);
393         
394   if(fSecondInputFileName!="") {
395           printf("Second Input File Name     =     %s\n", fSecondInputFileName.Data()) ;
396           printf("Second Input First Event   =     %d\n", fSecondInputFirstEvent) ;
397   }
398         
399   printf("Read Kine from, stack? %d, AOD ? %d \n", fReadStack, fReadAODMCParticles) ;
400   printf("Clean std AOD       =     %d\n", fCleanOutputStdAOD) ;
401   printf("Delta AOD File Name =     %s\n", fDeltaAODFileName.Data()) ;
402   printf("    \n") ;
403
404
405 //___________________________________________________
406 Bool_t AliCaloTrackReader::FillInputEvent(const Int_t iEntry, const char * currentFileName) {
407   //Fill the event counter and input lists that are needed, called by the analysis maker.
408
409   fEventNumber = iEntry;
410   fCurrentFileName = TString(currentFileName);
411   if(!fInputEvent) {
412           if(fDebug >= 0) printf("AliCaloTrackReader::FillInputEvent() - Input event not available, skip event analysis\n");
413           return kFALSE;
414   }
415   //Select events only fired by a certain trigger configuration if it is provided
416   Int_t eventType = 0;
417   if(fInputEvent->GetHeader())
418           eventType = ((AliVHeader*)fInputEvent->GetHeader())->GetEventType();
419   if( fFiredTriggerClassName  !="" && !fAnaLED){
420         if(eventType!=7)return kFALSE; //Only physics event, do not use for simulated events!!!
421     if(fDebug > 0) 
422       printf("AliCaloTrackReader::FillInputEvent() - FiredTriggerClass <%s>, selected class <%s>, compare name %d\n",
423              GetFiredTriggerClasses().Data(),fFiredTriggerClassName.Data(), GetFiredTriggerClasses().Contains(fFiredTriggerClassName));
424     if( !GetFiredTriggerClasses().Contains(fFiredTriggerClassName) ) return kFALSE;
425   }
426   else if(fAnaLED){
427 //        kStartOfRun =       1,    // START_OF_RUN
428 //        kEndOfRun =         2,    // END_OF_RUN
429 //        kStartOfRunFiles =  3,    // START_OF_RUN_FILES
430 //        kEndOfRunFiles =    4,    // END_OF_RUN_FILES
431 //        kStartOfBurst =     5,    // START_OF_BURST
432 //        kEndOfBurst =       6,    // END_OF_BURST
433 //        kPhysicsEvent =     7,    // PHYSICS_EVENT
434 //        kCalibrationEvent = 8,    // CALIBRATION_EVENT
435 //        kFormatError =      9,    // EVENT_FORMAT_ERROR
436 //        kStartOfData =      10,   // START_OF_DATA
437 //        kEndOfData =        11,   // END_OF_DATA
438 //        kSystemSoftwareTriggerEvent   = 12, // SYSTEM_SOFTWARE_TRIGGER_EVENT
439 //        kDetectorSoftwareTriggerEvent = 13  // DETECTOR_SOFTWARE_TRIGGER_EVENT
440          
441           if(eventType!=7 && fDebug > 1 )printf("AliCaloTrackReader::FillInputEvent() - DO LED, Event Type <%d>, 8 Calibration \n",  eventType);
442           if(eventType!=8)return kFALSE;
443   }
444                 
445   if(fOutputEvent && (fDataType != kAOD) && ((fOutputEvent->GetCaloClusters())->GetEntriesFast()!=0 ||(fOutputEvent->GetTracks())->GetEntriesFast()!=0)){
446     if (fFillCTS || fFillEMCAL || fFillPHOS) {
447       printf("AliCaloTrackReader::AODCaloClusters or AODTracks already filled by the filter, do not use the ESD reader, use the AOD reader, STOP\n");
448                   abort();
449     }
450   }
451
452   //In case of analysis of events with jets, skip those with jet pt > 5 pt hard 
453   if(fComparePtHardAndJetPt && GetStack()) {
454     if(!ComparePtHardAndJetPt()) return kFALSE ;
455   }
456
457   //In case of mixing events with other AOD file        
458   if(fDataType == kAOD && fSecondInputAODTree){
459          
460     if(fDebug > 1) 
461       printf("AliCaloTrackReader::FillInputEvent() - Get event %d from second input AOD file \n", iEntry+fSecondInputFirstEvent);
462     if(fSecondInputAODTree->GetEntriesFast() <= iEntry+fSecondInputFirstEvent) {
463       if(fSecondInputAODTree->GetEntriesFast() == iEntry+fSecondInputFirstEvent) 
464                          printf("AliCaloTrackReader::FillInputEvent() - Skip events from event %d, no more events in second AOD file \n", iEntry);
465       return kFALSE;
466     }
467     
468     //Get the Event
469     Int_t nbytes = fSecondInputAODTree->GetEvent(iEntry+fSecondInputFirstEvent);
470     if ( nbytes == 0 ) {//If nothing in AOD
471       printf("AliCaloTrackReader::FillInputEvent() - Nothing in Second AOD input, STOP\n");
472       abort() ; 
473     }
474     
475   }
476         
477   if(fFillEMCALCells) FillInputEMCALCells();
478   if(fFillPHOSCells)  FillInputPHOSCells();
479         
480   if(fFillCTS)   FillInputCTS();
481   if(fFillEMCAL) FillInputEMCAL();
482   if(fFillPHOS)  FillInputPHOS();
483
484         
485   return kTRUE ;
486 }
487
488 //__________________________________________________
489 void AliCaloTrackReader::ResetLists() {
490   //  Reset lists, called by the analysis maker 
491
492   if(fAODCTS)   fAODCTS -> Clear();
493   if(fAODEMCAL) fAODEMCAL -> Clear();
494   if(fAODPHOS)  fAODPHOS -> Clear();
495   if(fEMCALCells) fEMCALCells -> Clear();
496   if(fPHOSCells)  fPHOSCells -> Clear();
497   if(fCleanOutputStdAOD && fOutputEvent ){
498     //Only keep copied tracks and clusters if requested
499     fOutputEvent->GetTracks()      ->Clear();
500     fOutputEvent->GetCaloClusters()->Clear();
501   }
502
503 }