]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/CaloTrackCorrBase/AliCaloTrackReader.cxx
4ca5f67143e330b564e6f284f8c29d3d77351440
[u/mrichter/AliRoot.git] / PWG / CaloTrackCorrBase / 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 //-- Author: Gustavo Conesa (LNF-INFN) 
25 //////////////////////////////////////////////////////////////////////////////
26
27
28 // --- ROOT system ---
29 #include <TFile.h>
30 #include <TGeoManager.h>
31 #include <TGeoGlobalMagField.h>
32
33 // ---- ANALYSIS system ----
34 #include "AliMCEvent.h"
35 #include "AliAODMCHeader.h"
36 #include "AliGenPythiaEventHeader.h"
37 #include "AliESDEvent.h"
38 #include "AliAODEvent.h"
39 #include "AliVTrack.h"
40 #include "AliVParticle.h"
41 #include "AliMixedEvent.h"
42 #include "AliESDtrack.h"
43 #include "AliESDtrackCuts.h"
44 #include "AliTriggerAnalysis.h"
45 #include "AliESDVZERO.h"
46 #include "AliVCaloCells.h"
47 #include "AliMagF.h"
48
49 // ---- Detectors ----
50 #include "AliPHOSGeoUtils.h"
51 #include "AliEMCALGeometry.h"
52 #include "AliEMCALRecoUtils.h"
53
54 // ---- PartCorr ---
55 #include "AliCalorimeterUtils.h"
56 #include "AliCaloTrackReader.h"
57
58 ClassImp(AliCaloTrackReader)
59   
60   
61 //________________________________________
62 AliCaloTrackReader::AliCaloTrackReader() : 
63 TObject(),                   fEventNumber(-1), //fCurrentFileName(""),
64 fDataType(0),                fDebug(0), 
65 fFiducialCut(0x0),           fCheckFidCut(kFALSE), 
66 fComparePtHardAndJetPt(0),   fPtHardAndJetPtFactor(0),
67 fCTSPtMin(0),                fEMCALPtMin(0),                  fPHOSPtMin(0), 
68 fCTSPtMax(0),                fEMCALPtMax(0),                  fPHOSPtMax(0), 
69 fAODBranchList(0x0),
70 fCTSTracks(0x0),             fEMCALClusters(0x0),             fPHOSClusters(0x0),
71 fEMCALCells(0x0),            fPHOSCells(0x0),
72 fInputEvent(0x0),            fOutputEvent(0x0),fMC(0x0),
73 fFillCTS(0),                 fFillEMCAL(0),                   fFillPHOS(0),
74 fFillEMCALCells(0),          fFillPHOSCells(0), 
75 fRecalculateClusters(kFALSE),fSelectEmbeddedClusters(kFALSE),
76 fTrackStatus(0),             fTrackFilterMask(0),             fESDtrackCuts(0), 
77 fTrackMult(0),               fTrackMultEtaCut(0.8),
78 fReadStack(kFALSE),          fReadAODMCParticles(kFALSE), 
79 fDeltaAODFileName(""),       fFiredTriggerClassName(""),      fAnaLED(kFALSE),
80 fTaskName(""),               fCaloUtils(0x0), 
81 fMixedEvent(NULL),           fNMixedEvent(0),                 fVertex(NULL), 
82 fWriteOutputDeltaAOD(kFALSE),fOldAOD(kFALSE),                 fCaloFilterPatch(kFALSE),
83 fEMCALClustersListName(""),  fZvtxCut(0.),                    
84 fAcceptFastCluster(kFALSE),  fRemoveLEDEvents(kFALSE), 
85 fDoEventSelection(kFALSE),   fDoV0ANDEventSelection(kFALSE),  fUseEventsWithPrimaryVertex(kFALSE),
86 fTriggerAnalysis (0x0), 
87 fCentralityClass(""),        fCentralityOpt(0),
88 fEventPlaneMethod(""),       fImportGeometryFromFile(kFALSE), fImportGeometryFilePath("")
89
90 {
91   //Ctor
92
93   //Initialize parameters
94   InitParameters();
95 }
96
97 //_______________________________________
98 AliCaloTrackReader::~AliCaloTrackReader() 
99 {
100   //Dtor
101   
102   delete fFiducialCut ;
103         
104   if(fAODBranchList){
105     fAODBranchList->Delete();
106     delete fAODBranchList ;
107   }  
108   
109   if(fCTSTracks){
110     if(fDataType!=kMC)fCTSTracks->Clear() ; 
111     else              fCTSTracks->Delete() ; 
112     delete fCTSTracks ;
113   }
114   
115   if(fEMCALClusters){
116     if(fDataType!=kMC)fEMCALClusters->Clear("C") ; 
117     else              fEMCALClusters->Delete() ; 
118     delete fEMCALClusters ;
119   }
120   
121   if(fPHOSClusters){
122     if(fDataType!=kMC)fPHOSClusters->Clear("C") ; 
123     else              fPHOSClusters->Delete() ; 
124     delete fPHOSClusters ;
125   }
126   
127   if(fVertex){
128     for (Int_t i = 0; i < fNMixedEvent; i++) {
129       delete [] fVertex[i] ;
130       
131     }
132     delete [] fVertex ;
133         }
134   
135   delete fESDtrackCuts;
136   delete fTriggerAnalysis;
137   
138   //  Pointers not owned, done by the analysis frame
139   //  if(fInputEvent)  delete fInputEvent ;
140   //  if(fOutputEvent) delete fOutputEvent ;
141   //  if(fMC)          delete fMC ;  
142   //  Pointer not owned, deleted by maker
143   //  if (fCaloUtils) delete fCaloUtils ;
144   
145 }
146
147 //________________________________________________
148 Bool_t AliCaloTrackReader::ComparePtHardAndJetPt()
149 {
150   // Check the event, if the requested ptHard is much larger than the jet pT, then there is a problem.
151   // Only for PYTHIA.
152   if(!fReadStack) return kTRUE; //Information not filtered to AOD
153   
154   if(!strcmp(GetGenEventHeader()->ClassName(), "AliGenPythiaEventHeader")){
155     TParticle * jet =  0;
156     AliGenPythiaEventHeader* pygeh= (AliGenPythiaEventHeader*) GetGenEventHeader();
157     Int_t nTriggerJets =  pygeh->NTriggerJets();
158     Float_t ptHard = pygeh->GetPtHard();
159     
160     //if(fDebug > 1) printf("AliMCAnalysisUtils::PythiaEventHeader: Njets: %d, pT Hard %f\n",nTriggerJets, ptHard);
161     Float_t tmpjet[]={0,0,0,0};
162     for(Int_t ijet = 0; ijet< nTriggerJets; ijet++){
163       pygeh->TriggerJet(ijet, tmpjet);
164       jet = new TParticle(94, 21, -1, -1, -1, -1, tmpjet[0],tmpjet[1],tmpjet[2],tmpjet[3], 0,0,0,0);
165       //Compare jet pT and pt Hard
166       //if(fDebug > 1) printf("AliMCAnalysisUtils:: %d pycell jet pT %f\n",ijet, jet->Pt());
167       if(jet->Pt() > fPtHardAndJetPtFactor * ptHard) {
168         printf("AliMCAnalysisUtils::PythiaEventHeader: Njets: %d, pT Hard %2.2f, pycell jet pT %2.2f, rejection factor %1.1f\n",
169                nTriggerJets, ptHard, jet->Pt(), fPtHardAndJetPtFactor);
170         return kFALSE;
171       }
172     }
173     if(jet) delete jet; 
174   }
175   
176   return kTRUE ;
177   
178 }
179
180 //____________________________________________
181 AliStack* AliCaloTrackReader::GetStack() const 
182 {
183   //Return pointer to stack
184   if(fMC)
185     return fMC->Stack();
186   else{
187     if(fDebug > 1) printf("AliCaloTrackReader::GetStack() - Stack is not available\n"); 
188     return 0x0 ;
189   }
190 }
191
192 //__________________________________________________
193 TString AliCaloTrackReader::GetFiredTriggerClasses() 
194
195   // List of triggered classes in a TString
196
197   AliESDEvent* esdevent = dynamic_cast<AliESDEvent*> (GetInputEvent());
198   AliAODEvent* aodevent = dynamic_cast<AliAODEvent*> (GetInputEvent());
199     
200   if     (esdevent) return esdevent->GetFiredTriggerClasses();
201   else if(aodevent) return aodevent->GetFiredTriggerClasses();
202   else              return ""; // Mixed Event, MC event, does not have this trigger info
203
204 }
205
206 //______________________________________________
207 AliHeader* AliCaloTrackReader::GetHeader() const 
208 {
209   //Return pointer to header
210   if(fMC)
211     return fMC->Header();
212   else{
213     printf("AliCaloTrackReader::Header is not available\n"); 
214     return 0x0 ;
215   }
216 }
217
218 //______________________________________________________________
219 AliGenEventHeader* AliCaloTrackReader::GetGenEventHeader() const 
220 {
221   //Return pointer to Generated event header
222   if(fMC)
223     return fMC->GenEventHeader();
224   else{
225     printf("AliCaloTrackReader::GenEventHeader is not available\n"); 
226     return 0x0 ;
227   }
228 }
229
230 //____________________________________________________________________
231 TClonesArray* AliCaloTrackReader::GetAODMCParticles(Int_t input) const 
232 {
233   //Return list of particles in AOD. Do it for the corresponding input event.
234   
235   TClonesArray * rv = NULL ; 
236   if(fDataType == kAOD){
237     
238     if(input == 0){
239       //Normal input AOD
240       AliAODEvent * evt = dynamic_cast<AliAODEvent*> (fInputEvent) ;
241       if(evt)
242         rv = (TClonesArray*)evt->FindListObject("mcparticles");
243       else  
244         printf("AliCaloTrackReader::GetAODMCParticles() - wrong AOD input index? %d, or non existing tree? \n",input); 
245       
246     }  
247     
248   } else {
249     printf("AliCaloTrackReader::GetAODMCParticles() - Input are not AODs\n"); 
250   }
251   
252   return rv ; 
253 }
254
255 //___________________________________________________________________
256 AliAODMCHeader* AliCaloTrackReader::GetAODMCHeader(Int_t input) const 
257 {
258   //Return MC header in AOD. Do it for the corresponding input event.
259   AliAODMCHeader *mch = NULL;
260   if(fDataType == kAOD){
261     //Normal input AOD
262     if(input == 0) {
263       mch = (AliAODMCHeader*)((AliAODEvent*)fInputEvent)->FindListObject("mcheader");
264     }
265     //          //Second input AOD
266     //          else if(input == 1){ 
267     //       mch = (AliAODMCHeader*) fSecondInputAODEvent->FindListObject("mcheader");
268     //  }
269     else {
270       printf("AliCaloTrackReader::GetAODMCHeader() - wrong AOD input index, %d\n",input);
271     }
272   }
273   else {
274     printf("AliCaloTrackReader::GetAODMCHeader() - Input are not AODs\n");
275   }
276   
277   return mch;
278 }
279
280 //_____________________________
281 void AliCaloTrackReader::Init()
282 {
283   //Init reader. Method to be called in AliAnaPartCorrMaker
284
285   //printf(" AliCaloTrackReader::Init() %p \n",gGeoManager);
286
287   if(fReadStack && fReadAODMCParticles){
288     printf("AliCaloTrackReader::Init() - Cannot access stack and mcparticles at the same time, change them \n");
289     fReadStack          = kFALSE;
290     fReadAODMCParticles = kFALSE;
291   }
292   
293   // Init geometry, I do not like much to do it like this ...
294   if(fImportGeometryFromFile && !gGeoManager) {
295     printf("AliCaloTrackReader::Init() - Import geometry.root file\n");
296     TGeoManager::Import(Form("%s/geometry.root", fImportGeometryFilePath.Data())) ; // default need file "geometry.root" in local dir!!!!
297   }
298
299 }
300
301 //_______________________________________
302 void AliCaloTrackReader::InitParameters()
303 {
304   //Initialize the parameters of the analysis.
305   fDataType   = kESD ;
306   fCTSPtMin   = 0.1 ;
307   fEMCALPtMin = 0.1 ;
308   fPHOSPtMin  = 0.1 ;
309   fCTSPtMax   = 1000. ;
310   fEMCALPtMax = 1000. ;
311   fPHOSPtMax  = 1000. ;
312   
313   //Do not filter the detectors input by default.
314   fFillEMCAL      = kFALSE;
315   fFillPHOS       = kFALSE;
316   fFillCTS        = kFALSE;
317   fFillEMCALCells = kFALSE;
318   fFillPHOSCells  = kFALSE;
319   
320   fReadStack             = kFALSE; // Check in the constructor of the other readers if it was set or in the configuration file
321   fReadAODMCParticles    = kFALSE; // Check in the constructor of the other readers if it was set or in the configuration file
322   fDeltaAODFileName      = "deltaAODPartCorr.root";
323   fFiredTriggerClassName = "";
324   
325   fAcceptFastCluster = kTRUE;
326   fAnaLED            = kFALSE;
327   
328   //We want tracks fitted in the detectors:
329   //fTrackStatus=AliESDtrack::kTPCrefit;
330   //fTrackStatus|=AliESDtrack::kITSrefit; 
331   fTrackStatus     = 0;
332   fTrackFilterMask = 128; //For AODs, but what is the difference between fTrackStatus and fTrackFilterMask?
333   
334   fESDtrackCuts = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts(); //initialize with TPC only tracks 
335   
336   fV0ADC[0] = 0;   fV0ADC[1] = 0; 
337   fV0Mul[0] = 0;   fV0Mul[1] = 0; 
338   
339   fZvtxCut   = 10.;
340   
341   fNMixedEvent = 1;
342   
343   fPtHardAndJetPtFactor = 7;
344   
345   //Centrality
346   fCentralityClass  = "V0M";
347   fCentralityOpt    = 10;
348   fCentralityBin[0] = fCentralityBin[1]=-1;
349   
350   fEventPlaneMethod = "Q";
351
352   // Allocate memory (not sure this is the right place)
353   fCTSTracks       = new TObjArray();
354   fEMCALClusters   = new TObjArray();
355   fPHOSClusters    = new TObjArray(); 
356   fTriggerAnalysis = new AliTriggerAnalysis;
357   fAODBranchList   = new TList ;
358
359   fImportGeometryFromFile = kFALSE;
360   fImportGeometryFilePath = ".";
361   
362 }
363
364 //________________________________________________________
365 void AliCaloTrackReader::Print(const Option_t * opt) const
366 {
367   
368   //Print some relevant parameters set for the analysis
369   if(! opt)
370     return;
371   
372   printf("***** Print: %s %s ******\n", GetName(), GetTitle() ) ;
373   printf("Task name      : %s\n", fTaskName.Data()) ;
374   printf("Data type      : %d\n", fDataType) ;
375   printf("CTS Min pT     : %2.1f GeV/c\n", fCTSPtMin) ;
376   printf("EMCAL Min pT   : %2.1f GeV/c\n", fEMCALPtMin) ;
377   printf("PHOS Min pT    : %2.1f GeV/c\n", fPHOSPtMin) ;
378   printf("CTS Max pT     : %2.1f GeV/c\n", fCTSPtMax) ;
379   printf("EMCAL Max pT   : %2.1f GeV/c\n", fEMCALPtMax) ;
380   printf("PHOS Max pT    : %2.1f GeV/c\n", fPHOSPtMax) ;
381   printf("Use CTS         =     %d\n", fFillCTS) ;
382   printf("Use EMCAL       =     %d\n", fFillEMCAL) ;
383   printf("Use PHOS        =     %d\n", fFillPHOS) ;
384   printf("Use EMCAL Cells =     %d\n", fFillEMCALCells) ;
385   printf("Use PHOS  Cells =     %d\n", fFillPHOSCells) ;
386   printf("Track status    =     %d\n", (Int_t) fTrackStatus) ;
387   printf("Track filter mask (AODs) =  %d\n", (Int_t) fTrackFilterMask) ;
388   printf("Track Mult Eta Cut =  %d\n", (Int_t) fTrackMultEtaCut) ;
389   printf("Write delta AOD =     %d\n", fWriteOutputDeltaAOD) ;
390   printf("Recalculate Clusters = %d\n", fRecalculateClusters) ;
391   
392   if(fComparePtHardAndJetPt)
393           printf("Compare jet pt and pt hard to accept event, factor = %2.2f",fPtHardAndJetPtFactor);
394   
395   printf("Read Kine from, stack? %d, AOD ? %d \n", fReadStack, fReadAODMCParticles) ;
396   printf("Delta AOD File Name =     %s\n", fDeltaAODFileName.Data()) ;
397   printf("Centrality: Class %s, Option %d, Bin [%d,%d] \n", fCentralityClass.Data(),fCentralityOpt,fCentralityBin[0], fCentralityBin[1]) ;
398   
399   printf("    \n") ;
400   
401
402
403 //_________________________________________________________________________
404 Bool_t AliCaloTrackReader::FillInputEvent(const Int_t iEntry, 
405                                           const char * /*currentFileName*/) 
406 {
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   
420   if (GetFiredTriggerClasses().Contains("FAST")  && !GetFiredTriggerClasses().Contains("ALL") && !fAcceptFastCluster) {
421     if(fDebug > 0)  printf("AliCaloTrackReader::FillInputEvent - Do not count events from fast cluster, trigger name %s\n",fFiredTriggerClassName.Data());
422     return kFALSE;
423   }
424   
425   //-------------------------------------------------------------------------------------
426   // Reject event if large clusters with large energy
427   // Use only for LHC11a data for the moment, and if input is clusterizer V1 or V1+unfolding
428   // If clusterzer NxN or V2 it does not help
429   //-------------------------------------------------------------------------------------
430   if(fRemoveLEDEvents){
431     
432     //printf("Event %d\n",GetEventNumber());
433     for (Int_t i = 0; i < fInputEvent->GetNumberOfCaloClusters(); i++)
434     {
435       AliVCluster *clus = fInputEvent->GetCaloCluster(i);
436       if(clus->IsEMCAL()){               
437         if ((clus->E() > 500 && clus->GetNCells() > 200 ) || clus->GetNCells() > 200) {
438           Int_t absID = clus->GetCellsAbsId()[0];
439           Int_t sm = GetCaloUtils()->GetEMCALGeometry()->GetSuperModuleNumber(absID);
440           if(fDebug > 0) printf("AliCaloTrackReader::FillInputEvent - reject event %d with cluster : E %f, ncells %d, absId(0) %d, SM %d\n",GetEventNumber(),clus->E(),  clus->GetNCells(),absID, sm);
441           return kFALSE;
442         }
443       }
444     }
445     
446     // Count number of cells with energy larger than 0.1 in SM3, cut on this number
447     Int_t ncellsSM3 = 0;
448     Int_t ncellsSM4 = 0;
449     for(Int_t icell = 0; icell < fInputEvent->GetEMCALCells()->GetNumberOfCells(); icell++){
450       Int_t absID = fInputEvent->GetEMCALCells()->GetCellNumber(icell);
451       Int_t sm = GetCaloUtils()->GetEMCALGeometry()->GetSuperModuleNumber(absID);
452       if(fInputEvent->GetEMCALCells()->GetAmplitude(icell) > 0.1 && sm==3) ncellsSM3++;
453       if(fInputEvent->GetEMCALCells()->GetAmplitude(icell) > 0.1 && sm==4) ncellsSM4++;
454     }
455     
456     Int_t ncellcut = 21;
457     if(fFiredTriggerClassName.Contains("EMC")) ncellcut = 35;
458     
459     if(ncellsSM3 >= ncellcut || ncellsSM4 >= 100) {
460       if(fDebug > 0) printf(" AliCaloTrackReader::FillInputEvent() - reject event with ncells in SM3 %d and SM4 %d\n",ncellsSM3, ncellsSM4);
461       return kFALSE;
462     }
463   }// Remove LED events
464   
465   // Reject pure LED events?
466   if( fFiredTriggerClassName  !="" && !fAnaLED){
467     if(eventType!=7)
468       return kFALSE; //Only physics event, do not use for simulated events!!!
469     if(fDebug > 0) 
470       printf("AliCaloTrackReader::FillInputEvent() - FiredTriggerClass <%s>, selected class <%s>, compare name %d\n",
471              GetFiredTriggerClasses().Data(),fFiredTriggerClassName.Data(), GetFiredTriggerClasses().Contains(fFiredTriggerClassName));
472     if( !GetFiredTriggerClasses().Contains(fFiredTriggerClassName) ) return kFALSE;
473     else if(fDebug > 0) printf("AliCaloTrackReader::FillInputEvent() - Accepted triggered event\n");
474   }
475   else if(fAnaLED){
476     //    kStartOfRun =       1,    // START_OF_RUN
477     //    kEndOfRun =         2,    // END_OF_RUN
478     //    kStartOfRunFiles =  3,    // START_OF_RUN_FILES
479     //    kEndOfRunFiles =    4,    // END_OF_RUN_FILES
480     //    kStartOfBurst =     5,    // START_OF_BURST
481     //    kEndOfBurst =       6,    // END_OF_BURST
482     //    kPhysicsEvent =     7,    // PHYSICS_EVENT
483     //    kCalibrationEvent = 8,    // CALIBRATION_EVENT
484     //    kFormatError =      9,    // EVENT_FORMAT_ERROR
485     //    kStartOfData =      10,   // START_OF_DATA
486     //    kEndOfData =        11,   // END_OF_DATA
487     //    kSystemSoftwareTriggerEvent   = 12, // SYSTEM_SOFTWARE_TRIGGER_EVENT
488     //    kDetectorSoftwareTriggerEvent = 13  // DETECTOR_SOFTWARE_TRIGGER_EVENT
489     
490           if(eventType!=7 && fDebug > 1 )printf("AliCaloTrackReader::FillInputEvent() - DO LED, Event Type <%d>, 8 Calibration \n",  eventType);
491           if(eventType!=8)return kFALSE;
492   }
493   
494   //In case of analysis of events with jets, skip those with jet pt > 5 pt hard 
495   if(fComparePtHardAndJetPt && GetStack()) {
496     if(!ComparePtHardAndJetPt()) return kFALSE ;
497   }
498   
499   //Fill Vertex array
500   FillVertexArray();
501   //Reject events with Z vertex too large, only for SE analysis, if not, cut on the analysis code
502   if(!GetMixedEvent() && TMath::Abs(fVertex[0][2]) > fZvtxCut) return kFALSE;  
503   
504   //------------------------------------------------------
505   //Event rejection depending on vertex, pileup, v0and
506   //------------------------------------------------------
507   if(fDoEventSelection){
508     if(!fCaloFilterPatch){
509       //Do not analyze events with pileup
510       Bool_t bPileup = fInputEvent->IsPileupFromSPD(3, 0.8, 3., 2., 5.); //Default values, if not it does not compile
511       //Bool_t bPileup = event->IsPileupFromSPD(); 
512       if(bPileup) return kFALSE;
513       
514       if(fDoV0ANDEventSelection){
515         Bool_t bV0AND = kTRUE; 
516         AliESDEvent* esd = dynamic_cast<AliESDEvent*> (fInputEvent);
517         if(esd) 
518           bV0AND = fTriggerAnalysis->IsOfflineTriggerFired(esd, AliTriggerAnalysis::kV0AND);
519         //else bV0AND = //FIXME FOR AODs
520         if(!bV0AND) return kFALSE;
521       }
522       
523       if(fUseEventsWithPrimaryVertex && !CheckForPrimaryVertex()) return kFALSE;
524       
525     }//CaloFilter patch
526     else{ 
527       if(fInputEvent->GetNumberOfCaloClusters() > 0) {
528         AliVCluster * calo = fInputEvent->GetCaloCluster(0);
529         if(calo->GetNLabels() == 4){
530           Int_t * selection = calo->GetLabels();
531           Bool_t bPileup = selection[0];
532           if(bPileup) return kFALSE;
533           
534           Bool_t bGoodV = selection[1]; 
535           if(fUseEventsWithPrimaryVertex && !bGoodV) return kFALSE;
536           
537           if(fDoV0ANDEventSelection){
538             Bool_t bV0AND = selection[2]; 
539             if(!bV0AND) return kFALSE;
540           }
541           
542           fTrackMult = selection[3];
543           if(fTrackMult == 0) return kFALSE;
544         } else {
545           //First filtered AODs, track multiplicity stored there.  
546           fTrackMult = (Int_t) ((AliAODHeader*)fInputEvent->GetHeader())->GetCentrality();
547           if(fTrackMult == 0) return kFALSE;          
548         }
549       }//at least one cluster
550       else {
551         //printf("AliCaloTrackReader::FillInputEvent() - No clusters in event\n");
552         //Remove events with  vertex (0,0,0), bad vertex reconstruction
553         if(fUseEventsWithPrimaryVertex && TMath::Abs(fVertex[0][0]) < 1.e-6 && TMath::Abs(fVertex[0][1]) < 1.e-6 && TMath::Abs(fVertex[0][2]) < 1.e-6) return kFALSE;
554         
555         //First filtered AODs, track multiplicity stored there.  
556         fTrackMult = (Int_t) ((AliAODHeader*)fInputEvent->GetHeader())->GetCentrality();
557         if(fTrackMult == 0) return kFALSE;
558       }// no cluster
559     }// CaloFileter patch
560   }// Event selection
561   //------------------------------------------------------
562   
563   //Check if there is a centrality value, PbPb analysis, and if a centrality bin selection is requested
564   //If we need a centrality bin, we select only those events in the corresponding bin.
565   if(GetCentrality() && fCentralityBin[0]>=0 && fCentralityBin[1]>=0 && fCentralityOpt==100){
566     Int_t cen = GetEventCentrality();
567     if(cen > fCentralityBin[1] || cen < fCentralityBin[0]) return kFALSE; //reject events out of bin.
568   }
569   
570   //Fill the arrays with cluster/tracks/cells data
571   if(fFillEMCALCells) 
572     FillInputEMCALCells();
573   if(fFillPHOSCells)  
574     FillInputPHOSCells();
575         
576   if(fFillCTS){   
577     FillInputCTS();
578     //Accept events with at least one track
579     if(fTrackMult == 0 && fDoEventSelection) return kFALSE;
580   }
581   
582   if(fFillEMCAL) 
583     FillInputEMCAL();
584   if(fFillPHOS)  
585     FillInputPHOS();
586   
587   FillInputVZERO();
588         
589   return kTRUE ;
590 }
591
592 //___________________________________
593 void AliCaloTrackReader::ResetLists() 
594 {
595   //  Reset lists, called by the analysis maker 
596   
597   if(fCTSTracks)       fCTSTracks     -> Clear();
598   if(fEMCALClusters)   fEMCALClusters -> Clear("C");
599   if(fPHOSClusters)    fPHOSClusters  -> Clear("C");
600   
601   fV0ADC[0] = 0;   fV0ADC[1] = 0; 
602   fV0Mul[0] = 0;   fV0Mul[1] = 0; 
603   
604 }
605
606 //____________________________________________________________
607 void AliCaloTrackReader::SetInputEvent(AliVEvent* const input)  
608 {
609   fInputEvent  = input;
610   fMixedEvent = dynamic_cast<AliMixedEvent*>(GetInputEvent()) ; 
611   if (fMixedEvent) {
612     fNMixedEvent = fMixedEvent->GetNumberOfEvents() ; 
613   }
614   
615   //Delete previous vertex
616   if(fVertex){
617     for (Int_t i = 0; i < fNMixedEvent; i++) {
618       delete [] fVertex[i] ; 
619     }
620     delete [] fVertex ;
621   }
622   
623   fVertex = new Double_t*[fNMixedEvent] ; 
624   for (Int_t i = 0; i < fNMixedEvent; i++) {
625     fVertex[i] = new Double_t[3] ; 
626     fVertex[i][0] = 0.0 ; 
627     fVertex[i][1] = 0.0 ; 
628     fVertex[i][2] = 0.0 ; 
629   }
630 }
631
632 //__________________________________________________
633 Int_t AliCaloTrackReader::GetEventCentrality() const 
634 {
635   //Return current event centrality
636   
637   if(GetCentrality()){
638     if(fCentralityOpt==100)     return (Int_t) GetCentrality()->GetCentralityPercentile(fCentralityClass); // 100 bins max
639     else if(fCentralityOpt==10) return GetCentrality()->GetCentralityClass10(fCentralityClass);// 10 bins max
640     else if(fCentralityOpt==20) return GetCentrality()->GetCentralityClass5(fCentralityClass); // 20 bins max
641     else {
642       printf("AliAnaPartCorrBaseClass::Unknown centrality option %d, use 10, 20 or 100\n",fCentralityOpt);
643       return 0;
644     } 
645   }
646   else return 0;
647   
648 }
649
650 //__________________________________________________________
651 void AliCaloTrackReader::GetVertex(Double_t vertex[3]) const 
652 {
653   //Return vertex position to be used for single event analysis
654   vertex[0]=fVertex[0][0];  
655   vertex[1]=fVertex[0][1];  
656   vertex[2]=fVertex[0][2];
657 }
658
659 //____________________________________________________________
660 void AliCaloTrackReader::GetVertex(Double_t vertex[3], 
661                                    const Int_t evtIndex) const 
662 {
663   //Return vertex position for mixed event, recover the vertex in a particular event.
664   
665   vertex[0]=fVertex[evtIndex][0];  vertex[1]=fVertex[evtIndex][1];  vertex[2]=fVertex[evtIndex][2];
666   
667 }
668
669 //________________________________________
670 void AliCaloTrackReader::FillVertexArray() 
671 {
672   
673   //Fill data member with vertex
674   //In case of Mixed event, multiple vertices
675   
676   //Delete previous vertex
677   if(fVertex){
678     for (Int_t i = 0; i < fNMixedEvent; i++) {
679       delete [] fVertex[i] ; 
680     }
681     delete [] fVertex ;  
682   }
683   
684   fVertex = new Double_t*[fNMixedEvent] ; 
685   for (Int_t i = 0; i < fNMixedEvent; i++) {
686     fVertex[i] = new Double_t[3] ; 
687     fVertex[i][0] = 0.0 ; 
688     fVertex[i][1] = 0.0 ; 
689     fVertex[i][2] = 0.0 ; 
690   }          
691   
692   if (!fMixedEvent) { //Single event analysis
693     if(fDataType!=kMC){
694       
695       if(fInputEvent->GetPrimaryVertex()){
696         fInputEvent->GetPrimaryVertex()->GetXYZ(fVertex[0]); 
697       }
698       else {
699         printf("AliCaloTrackReader::FillVertexArray() - NULL primary vertex\n");
700         fVertex[0][0]=0.;   fVertex[0][1]=0.;   fVertex[0][2]=0.;
701       }//Primary vertex pointer do not exist
702       
703     } else {//MC read event 
704       fVertex[0][0]=0.;   fVertex[0][1]=0.;   fVertex[0][2]=0.;
705     }
706     
707     if(fDebug > 1)
708       printf("AliCaloTrackReader::FillVertexArray() - Single Event Vertex : %f,%f,%f\n",fVertex[0][0],fVertex[0][1],fVertex[0][2]);
709     
710   } else { // MultiEvent analysis
711     for (Int_t iev = 0; iev < fNMixedEvent; iev++) {
712       if (fMixedEvent->GetVertexOfEvent(iev))
713         fMixedEvent->GetVertexOfEvent(iev)->GetXYZ(fVertex[iev]);
714       else { // no vertex found !!!!
715         AliWarning("No vertex found");
716       }
717       
718       if(fDebug > 1)
719         printf("AliCaloTrackReader::FillVertexArray() - Multi Event %d Vertex : %f,%f,%f\n",iev,fVertex[iev][0],fVertex[iev][1],fVertex[iev][2]);
720       
721     }
722   }
723   
724 }
725
726 //_____________________________________
727 void AliCaloTrackReader::FillInputCTS() 
728 {
729   //Return array with Central Tracking System (CTS) tracks
730   
731   if(fDebug > 2 ) printf("AliCaloTrackReader::FillInputCTS()\n");
732   
733   Int_t nTracks   = fInputEvent->GetNumberOfTracks() ;
734   Double_t p[3];
735   fTrackMult = 0;
736   Int_t nstatus = 0;
737   for (Int_t itrack =  0; itrack <  nTracks; itrack++) {////////////// track loop
738     AliVTrack * track = (AliVTrack*)fInputEvent->GetTrack(itrack) ; // retrieve track from esd
739     
740     //Select tracks under certain conditions, TPCrefit, ITSrefit ... check the set bits
741     if (fTrackStatus && !((track->GetStatus() & fTrackStatus) == fTrackStatus)) 
742       continue ;
743     
744     nstatus++;
745     
746     if     (fDataType==kESD && !fESDtrackCuts->AcceptTrack((AliESDtrack*)track))
747     {
748       continue;
749     }
750     else if(fDataType==kAOD)
751     {
752       AliAODTrack *aodtrack = dynamic_cast <AliAODTrack*>(track);
753       if(aodtrack){
754         if(fDebug > 2 ) 
755           printf("AliCaloTrackReader::FillInputCTS():AOD track type: %c \n", aodtrack->GetType());
756         if (fDataType!=kMC && aodtrack->TestFilterBit(fTrackFilterMask)==kFALSE) continue;
757         if(aodtrack->GetType()!=AliAODTrack::kPrimary) continue;
758       }
759     }
760     
761     //Count the tracks in eta < 0.9
762     //printf("Eta %f cut  %f\n",TMath::Abs(track->Eta()),fTrackMultEtaCut);
763     if(TMath::Abs(track->Eta())< fTrackMultEtaCut) fTrackMult++;
764     
765     track->GetPxPyPz(p) ;
766     TLorentzVector momentum(p[0],p[1],p[2],0);
767     
768     if(fCTSPtMin < momentum.Pt() && fCTSPtMax > momentum.Pt()){
769       
770       if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"CTS")) 
771         continue;
772       
773       if(fDebug > 2 && momentum.Pt() > 0.1) 
774         printf("AliCaloTrackReader::FillInputCTS() - Selected tracks E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
775                momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
776       
777       if (fMixedEvent) {
778         track->SetID(itrack);
779       }
780       
781       fCTSTracks->Add(track);        
782       
783     }//Pt and Fiducial cut passed. 
784   }// track loop
785         
786   //fCTSTracksNormalInputEntries = fCTSTracks->GetEntriesFast();
787   if(fDebug > 1) 
788     printf("AliCaloTrackReader::FillInputCTS()   - aod entries %d, input tracks %d, pass status %d, multipliticy %d\n", fCTSTracks->GetEntriesFast(), nTracks, nstatus, fTrackMult);//fCTSTracksNormalInputEntries);
789   
790 }
791
792 //__________________________________________________________________
793 void AliCaloTrackReader::FillInputEMCALAlgorithm(AliVCluster * clus, 
794                                                  const Int_t iclus) 
795 {
796   //Fill the EMCAL data in the array, do it
797   
798   Int_t vindex = 0 ;  
799   if (fMixedEvent) 
800     vindex = fMixedEvent->EventIndexForCaloCluster(iclus);
801   
802   //Reject clusters with bad channels, close to borders and exotic;
803   if(!GetCaloUtils()->GetEMCALRecoUtils()->IsGoodCluster(clus,GetCaloUtils()->GetEMCALGeometry(),GetEMCALCells(),fInputEvent->GetBunchCrossNumber())) return;
804   
805   //Mask all cells in collumns facing ALICE thick material if requested
806   if(GetCaloUtils()->GetNMaskCellColumns()){
807     Int_t absId   = -1;
808     Int_t iSupMod = -1;
809     Int_t iphi    = -1;
810     Int_t ieta    = -1;
811     Bool_t shared = kFALSE;
812     GetCaloUtils()->GetEMCALRecoUtils()->GetMaxEnergyCell(GetCaloUtils()->GetEMCALGeometry(), GetEMCALCells(),clus,absId,iSupMod,ieta,iphi,shared);
813     if(GetCaloUtils()->MaskFrameCluster(iSupMod, ieta)) return;
814   }
815   
816   if(fSelectEmbeddedClusters){
817     if(clus->GetNLabels()==0 || clus->GetLabel() < 0) return;
818     //else printf("Embedded cluster,  %d, n label %d label %d  \n",iclus,clus->GetNLabels(),clus->GetLabel());
819   }
820   
821   //Float_t pos[3];
822   //clus->GetPosition(pos);
823   //printf("Before Corrections: e %f, x %f, y %f, z %f\n",clus->E(),pos[0],pos[1],pos[2]);
824   
825   if(fRecalculateClusters){
826     //Recalibrate the cluster energy 
827     if(GetCaloUtils()->IsRecalibrationOn()) {
828       
829       Float_t energy = GetCaloUtils()->RecalibrateClusterEnergy(clus, GetEMCALCells());
830       
831       clus->SetE(energy);
832       //printf("Recalibrated Energy %f\n",clus->E());  
833       
834       GetCaloUtils()->RecalculateClusterShowerShapeParameters(GetEMCALCells(),clus);
835       GetCaloUtils()->RecalculateClusterPID(clus);
836       
837     } // recalculate E
838     
839     //Recalculate distance to bad channels, if new list of bad channels provided
840     GetCaloUtils()->RecalculateClusterDistanceToBadChannel(GetEMCALCells(),clus);
841     
842     //Recalculate cluster position
843     if(GetCaloUtils()->IsRecalculationOfClusterPositionOn()){
844       GetCaloUtils()->RecalculateClusterPosition(GetEMCALCells(),clus); 
845       //clus->GetPosition(pos);
846       //printf("After  Corrections: e %f, x %f, y %f, z %f\n",clus->E(),pos[0],pos[1],pos[2]);
847     }
848   }
849   
850   // Recalculate TOF
851   if(GetCaloUtils()->GetEMCALRecoUtils()->IsTimeRecalibrationOn()) {
852     
853     Double_t tof      = clus->GetTOF();
854     Float_t  frac     =-1;
855     Int_t    absIdMax = GetCaloUtils()->GetMaxEnergyCell(fEMCALCells, clus,frac);
856     
857     if(fDataType==AliCaloTrackReader::kESD){ 
858       tof = fEMCALCells->GetCellTime(absIdMax);
859     }
860     
861     GetCaloUtils()->GetEMCALRecoUtils()->RecalibrateCellTime(absIdMax,fInputEvent->GetBunchCrossNumber(),tof);
862     
863     clus->SetTOF(tof);
864     
865   }// Time recalibration    
866   
867   //Correct non linearity
868   if(GetCaloUtils()->IsCorrectionOfClusterEnergyOn()){
869     GetCaloUtils()->CorrectClusterEnergy(clus) ;
870     //printf("Linearity Corrected Energy %f\n",clus->E());  
871     
872     //In case of MC analysis, to match resolution/calibration in real data
873     Float_t rdmEnergy = GetCaloUtils()->GetEMCALRecoUtils()->SmearClusterEnergy(clus);
874     // printf("\t Energy %f, smeared %f\n", clus->E(),rdmEnergy);
875     clus->SetE(rdmEnergy);
876   }
877   
878   TLorentzVector momentum ;
879   
880   clus->GetMomentum(momentum, fVertex[vindex]);      
881   
882   if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"EMCAL")) return;
883   
884   if(fEMCALPtMin > momentum.E() || fEMCALPtMax < momentum.E())         return;
885   
886   if(fDebug > 2 && momentum.E() > 0.1) 
887     printf("AliCaloTrackReader::FillInputEMCAL() - Selected clusters E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
888            momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
889   
890   if (fMixedEvent) 
891     clus->SetID(iclus) ; 
892   
893   fEMCALClusters->Add(clus);    
894   
895 }
896
897 //_______________________________________
898 void AliCaloTrackReader::FillInputEMCAL() 
899 {
900   //Return array with EMCAL clusters in aod format
901   
902   if(fDebug > 2 ) printf("AliCaloTrackReader::FillInputEMCAL()\n");
903   
904   // First recalibrate cells, time or energy
905   //  if(GetCaloUtils()->IsRecalibrationOn())
906   //    GetCaloUtils()->GetEMCALRecoUtils()->RecalibrateCells(GetCaloUtils()->GetEMCALGeometry(), 
907   //                                                          GetEMCALCells(), 
908   //                                                          fInputEvent->GetBunchCrossNumber());
909   
910   //Loop to select clusters in fiducial cut and fill container with aodClusters
911   if(fEMCALClustersListName==""){
912     Int_t nclusters = fInputEvent->GetNumberOfCaloClusters();
913     for (Int_t iclus =  0; iclus <  nclusters; iclus++) {
914       AliVCluster * clus = 0;
915       if ( (clus = fInputEvent->GetCaloCluster(iclus)) ) {
916         if (IsEMCALCluster(clus)){          
917           FillInputEMCALAlgorithm(clus, iclus);
918         }//EMCAL cluster
919       }// cluster exists
920     }// cluster loop
921     
922     //Recalculate track matching
923     GetCaloUtils()->RecalculateClusterTrackMatching(fInputEvent);
924     
925   }//Get the clusters from the input event
926   else {
927     TClonesArray * clusterList = 0x0; 
928     if(fOutputEvent) clusterList = dynamic_cast<TClonesArray*> (fOutputEvent->FindListObject(fEMCALClustersListName));
929     if(!clusterList){
930       //printf("AliCaloTrackReader::FillInputEMCAL() - Wrong name of list with clusters? Try input event <%s>\n",fEMCALClustersListName.Data());
931       //List not in output event, try input event
932       clusterList = dynamic_cast<TClonesArray*> (fInputEvent->FindListObject(fEMCALClustersListName));
933       if(!clusterList){
934         printf("AliCaloTrackReader::FillInputEMCAL() - Wrong name of list with clusters?  <%s>\n",fEMCALClustersListName.Data());
935         return;
936       }
937     }
938     
939     Int_t nclusters = clusterList->GetEntriesFast();
940     for (Int_t iclus =  0; iclus <  nclusters; iclus++) {
941       AliVCluster * clus = dynamic_cast<AliVCluster*> (clusterList->At(iclus));
942       //printf("E %f\n",clus->E());
943       if (clus) FillInputEMCALAlgorithm(clus, iclus);
944       else printf("AliCaloTrackReader::FillInputEMCAL() - Null cluster in list!\n");
945       
946     }// cluster loop
947     
948     // Recalculate track matching, not necessary if already done in the reclusterization task.
949     // in case it was not done ...
950     GetCaloUtils()->RecalculateClusterTrackMatching(fInputEvent,clusterList);
951     
952   }
953   
954   if(fDebug > 1) printf("AliCaloTrackReader::FillInputEMCAL() - aod entries %d\n",  fEMCALClusters->GetEntriesFast());
955   
956 }
957
958 //______________________________________
959 void AliCaloTrackReader::FillInputPHOS() 
960 {
961   //Return array with PHOS clusters in aod format
962   
963   if(fDebug > 2 ) printf("AliCaloTrackReader::FillInputPHOS()\n");
964   
965   //Loop to select clusters in fiducial cut and fill container with aodClusters
966   Int_t nclusters = fInputEvent->GetNumberOfCaloClusters();
967   for (Int_t iclus = 0; iclus < nclusters; iclus++) {
968     AliVCluster * clus = 0;
969     if ( (clus = fInputEvent->GetCaloCluster(iclus)) ) {
970       if (IsPHOSCluster(clus)){
971         //Check if the cluster contains any bad channel and if close to calorimeter borders
972         Int_t vindex = 0 ;  
973         if (fMixedEvent) 
974           vindex = fMixedEvent->EventIndexForCaloCluster(iclus);
975         if( GetCaloUtils()->ClusterContainsBadChannel("PHOS",clus->GetCellsAbsId(), clus->GetNCells())) 
976           continue;
977         if(!GetCaloUtils()->CheckCellFiducialRegion(clus, fInputEvent->GetPHOSCells(), fInputEvent, vindex)) 
978           continue;
979         
980         if(fRecalculateClusters){
981           
982           //Recalibrate the cluster energy 
983           if(GetCaloUtils()->IsRecalibrationOn()) {
984             Float_t energy = GetCaloUtils()->RecalibrateClusterEnergy(clus, (AliAODCaloCells*)GetPHOSCells());
985             clus->SetE(energy);
986           }
987           
988         }
989         
990         TLorentzVector momentum ;
991         
992         clus->GetMomentum(momentum, fVertex[vindex]);      
993         
994         if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"PHOS")) continue;
995         
996         if(fPHOSPtMin > momentum.E() || fPHOSPtMax < momentum.E())          continue;
997         
998         if(fDebug > 2 && momentum.E() > 0.1) 
999           printf("AliCaloTrackReader::FillInputPHOS() - Selected clusters E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
1000                  momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());        
1001         
1002         
1003         if (fMixedEvent) {
1004           clus->SetID(iclus) ; 
1005         }              
1006         
1007         fPHOSClusters->Add(clus);       
1008         
1009       }//PHOS cluster
1010     }//cluster exists
1011   }//esd cluster loop
1012   
1013   if(fDebug > 1) printf("AliCaloTrackReader::FillInputPHOS()  - aod entries %d\n",  fPHOSClusters->GetEntriesFast());
1014   
1015 }
1016
1017 //____________________________________________
1018 void AliCaloTrackReader::FillInputEMCALCells() 
1019 {
1020   //Return array with EMCAL cells in aod format
1021   
1022   fEMCALCells = fInputEvent->GetEMCALCells(); 
1023   
1024 }
1025
1026 //___________________________________________
1027 void AliCaloTrackReader::FillInputPHOSCells() 
1028 {
1029   //Return array with PHOS cells in aod format
1030   
1031   fPHOSCells = fInputEvent->GetPHOSCells(); 
1032   
1033 }
1034
1035 //_______________________________________
1036 void AliCaloTrackReader::FillInputVZERO()
1037 {
1038   //Fill VZERO information in data member, add all the channels information.
1039   AliVVZERO* v0 = fInputEvent->GetVZEROData();
1040   //printf("Init V0: ADC (%d,%d), Multiplicity (%d,%d) \n",fV0ADC[0],fV0ADC[1],fV0Mul[0],fV0Mul[1]);
1041   
1042   if (v0) 
1043   {
1044     AliESDVZERO* esdV0 = dynamic_cast<AliESDVZERO*> (v0);
1045     for (Int_t i = 0; i < 32; i++)
1046     {
1047       if(esdV0){//Only available in ESDs
1048         fV0ADC[0] += (Int_t)esdV0->GetAdcV0C(i);
1049         fV0ADC[1] += (Int_t)esdV0->GetAdcV0A(i);
1050       }
1051       fV0Mul[0] += (Int_t)v0->GetMultiplicityV0C(i);
1052       fV0Mul[1] += (Int_t)v0->GetMultiplicityV0A(i);
1053     }
1054     if(fDebug > 0)
1055       printf("V0: ADC (%d,%d), Multiplicity (%d,%d) \n",fV0ADC[0],fV0ADC[1],fV0Mul[0],fV0Mul[1]);
1056   }
1057   else
1058   {
1059     if(fDebug > 0)
1060       printf("Cannot retrieve V0 ESD! Run w/ null V0 charges\n ");
1061   }
1062 }
1063
1064
1065 //___________________________________________________________________
1066 Bool_t AliCaloTrackReader::IsEMCALCluster(AliVCluster* cluster) const 
1067 {
1068   // Check if it is a cluster from EMCAL. For old AODs cluster type has
1069   // different number and need to patch here
1070   
1071   if(fDataType==kAOD && fOldAOD)
1072   {
1073     if (cluster->GetType() == 2) return kTRUE;
1074     else                         return kFALSE;
1075   }
1076   else 
1077   {
1078     return cluster->IsEMCAL();
1079   }
1080   
1081 }
1082
1083 //___________________________________________________________________
1084 Bool_t AliCaloTrackReader::IsPHOSCluster(AliVCluster * cluster) const 
1085 {
1086   //Check if it is a cluster from PHOS.For old AODs cluster type has
1087   // different number and need to patch here
1088   
1089   if(fDataType==kAOD && fOldAOD)
1090   {
1091     Int_t type = cluster->GetType();
1092     if (type == 0 || type == 1) return kTRUE;
1093     else                        return kFALSE;
1094   }
1095   else 
1096   {
1097     return cluster->IsPHOS();
1098   }
1099   
1100 }
1101
1102 //________________________________________________
1103 Bool_t AliCaloTrackReader::CheckForPrimaryVertex()
1104 {
1105   //Check if the vertex was well reconstructed, copy from V0Reader of conversion group
1106   //Only for ESDs ...
1107   
1108   AliESDEvent * event = dynamic_cast<AliESDEvent*> (fInputEvent);
1109   if(!event) return kTRUE;
1110   
1111   if(event->GetPrimaryVertexTracks()->GetNContributors() > 0) {
1112     return kTRUE;
1113   }
1114   
1115   if(event->GetPrimaryVertexTracks()->GetNContributors() < 1) {
1116     // SPD vertex
1117     if(event->GetPrimaryVertexSPD()->GetNContributors() > 0) {
1118       //cout<<"spd vertex type::"<< fESDEvent->GetPrimaryVertex()->GetName() << endl;
1119       return kTRUE;
1120       
1121     }
1122     if(event->GetPrimaryVertexSPD()->GetNContributors() < 1) {
1123       //      cout<<"bad vertex type::"<< event->GetPrimaryVertex()->GetName() << endl;
1124       return kFALSE;
1125     }
1126   }
1127   
1128   return kFALSE;  
1129   
1130 }
1131
1132 //____________________________________________________________
1133 void  AliCaloTrackReader::SetTrackCuts(AliESDtrackCuts * cuts) 
1134
1135   // Set Track cuts
1136   
1137   if(fESDtrackCuts) delete fESDtrackCuts ;
1138   
1139   fESDtrackCuts = cuts ; 
1140   
1141 }                 
1142