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