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