]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/PartCorrBase/AliCaloTrackReader.cxx
add protection in reader against missing output event when the non standard cluster...
[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), fCTSPtMax(1000), fEMCALPtMax(1000),fPHOSPtMax(1000), 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     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), fUseEventsWithPrimaryVertex(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.1 ;
312   fEMCALPtMin = 0.1 ;
313   fPHOSPtMin  = 0.1 ;
314   fCTSPtMax   = 1000. ;
315   fEMCALPtMax = 1000. ;
316   fPHOSPtMax  = 1000. ;
317   
318   //Do not filter the detectors input by default.
319   fFillEMCAL      = kFALSE;
320   fFillPHOS       = kFALSE;
321   fFillCTS        = kFALSE;
322   fFillEMCALCells = kFALSE;
323   fFillPHOSCells  = kFALSE;
324
325   //fSecondInputFileName   = "" ;
326   //fSecondInputFirstEvent = 0 ;
327   fReadStack             = kFALSE; // Check in the constructor of the other readers if it was set or in the configuration file
328   fReadAODMCParticles    = kFALSE; // Check in the constructor of the other readers if it was set or in the configuration file
329   fDeltaAODFileName      = "deltaAODPartCorr.root";
330   fFiredTriggerClassName      = "";
331                 
332   fAnaLED = kFALSE;
333   
334   //We want tracks fitted in the detectors:
335   //fTrackStatus=AliESDtrack::kTPCrefit;
336   //fTrackStatus|=AliESDtrack::kITSrefit;  
337   
338   fESDtrackCuts = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010();
339
340   fV0ADC[0] = 0;   fV0ADC[1] = 0; 
341   fV0Mul[0] = 0;   fV0Mul[1] = 0; 
342
343   fZvtxCut   = 10.;
344   
345   //Centrality
346   fCentralityBin[0]=fCentralityBin[1]=-1;
347   
348   //Cluster smearing
349   fSmearClusterEnergy   = kFALSE;
350   fSmearClusterParam[0] = 0.07; // * sqrt E term
351   fSmearClusterParam[1] = 0.02; // * E term
352   fSmearClusterParam[2] = 0.00; // constant
353   
354 }
355
356 //________________________________________________________________
357 void AliCaloTrackReader::Print(const Option_t * opt) const
358 {
359
360   //Print some relevant parameters set for the analysis
361   if(! opt)
362     return;
363
364   printf("***** Print: %s %s ******\n", GetName(), GetTitle() ) ;
365   printf("Task name      : %s\n", fTaskName.Data()) ;
366   printf("Data type      : %d\n", fDataType) ;
367   printf("CTS Min pT     : %2.1f GeV/c\n", fCTSPtMin) ;
368   printf("EMCAL Min pT   : %2.1f GeV/c\n", fEMCALPtMin) ;
369   printf("PHOS Min pT    : %2.1f GeV/c\n", fPHOSPtMin) ;
370   printf("CTS Max pT     : %2.1f GeV/c\n", fCTSPtMax) ;
371   printf("EMCAL Max pT   : %2.1f GeV/c\n", fEMCALPtMax) ;
372   printf("PHOS Max pT    : %2.1f GeV/c\n", fPHOSPtMax) ;
373   printf("Use CTS         =     %d\n", fFillCTS) ;
374   printf("Use EMCAL       =     %d\n", fFillEMCAL) ;
375   printf("Use PHOS        =     %d\n", fFillPHOS) ;
376   printf("Use EMCAL Cells =     %d\n", fFillEMCALCells) ;
377   printf("Use PHOS  Cells =     %d\n", fFillPHOSCells) ;
378   printf("Track status    =     %d\n", (Int_t) fTrackStatus) ;
379   printf("Track Mult Eta Cut =  %d\n", (Int_t) fTrackMultEtaCut) ;
380   printf("Write delta AOD =     %d\n", fWriteOutputDeltaAOD) ;
381
382   if(fComparePtHardAndJetPt)
383           printf("Compare jet pt and pt hard to accept event, factor = %2.2f",fPtHardAndJetPtFactor);
384         
385 //  if(fSecondInputFileName!="") {
386 //        printf("Second Input File Name     =     %s\n", fSecondInputFileName.Data()) ;
387 //        printf("Second Input First Event   =     %d\n", fSecondInputFirstEvent) ;
388 //  }
389         
390   printf("Read Kine from, stack? %d, AOD ? %d \n", fReadStack, fReadAODMCParticles) ;
391   printf("Delta AOD File Name =     %s\n", fDeltaAODFileName.Data()) ;
392   printf("Centrality: Class %s, Option %d, Bin [%d,%d] \n", fCentralityClass.Data(),fCentralityOpt,fCentralityBin[0], fCentralityBin[1]) ;
393
394   printf("    \n") ;
395   
396
397
398 //___________________________________________________
399 Bool_t AliCaloTrackReader::FillInputEvent(const Int_t iEntry, const char * /*currentFileName*/) {
400   //Fill the event counter and input lists that are needed, called by the analysis maker.
401
402   fEventNumber = iEntry;
403   //fCurrentFileName = TString(currentFileName);
404   if(!fInputEvent) {
405           if(fDebug >= 0) printf("AliCaloTrackReader::FillInputEvent() - Input event not available, skip event analysis\n");
406           return kFALSE;
407   }
408   //Select events only fired by a certain trigger configuration if it is provided
409   Int_t eventType = 0;
410   if(fInputEvent->GetHeader())
411           eventType = ((AliVHeader*)fInputEvent->GetHeader())->GetEventType();
412
413   if( fFiredTriggerClassName  !="" && !fAnaLED){
414     if(eventType!=7)
415       return kFALSE; //Only physics event, do not use for simulated events!!!
416     if(fDebug > 0) 
417       printf("AliCaloTrackReader::FillInputEvent() - FiredTriggerClass <%s>, selected class <%s>, compare name %d\n",
418              GetFiredTriggerClasses().Data(),fFiredTriggerClassName.Data(), GetFiredTriggerClasses().Contains(fFiredTriggerClassName));
419     if( !GetFiredTriggerClasses().Contains(fFiredTriggerClassName) ) return kFALSE;
420     else if(fDebug > 0) printf("AliCaloTrackReader::FillInputEvent() - Accepted triggered event\n");
421   }
422   else if(fAnaLED){
423 //        kStartOfRun =       1,    // START_OF_RUN
424 //        kEndOfRun =         2,    // END_OF_RUN
425 //        kStartOfRunFiles =  3,    // START_OF_RUN_FILES
426 //        kEndOfRunFiles =    4,    // END_OF_RUN_FILES
427 //        kStartOfBurst =     5,    // START_OF_BURST
428 //        kEndOfBurst =       6,    // END_OF_BURST
429 //        kPhysicsEvent =     7,    // PHYSICS_EVENT
430 //        kCalibrationEvent = 8,    // CALIBRATION_EVENT
431 //        kFormatError =      9,    // EVENT_FORMAT_ERROR
432 //        kStartOfData =      10,   // START_OF_DATA
433 //        kEndOfData =        11,   // END_OF_DATA
434 //        kSystemSoftwareTriggerEvent   = 12, // SYSTEM_SOFTWARE_TRIGGER_EVENT
435 //        kDetectorSoftwareTriggerEvent = 13  // DETECTOR_SOFTWARE_TRIGGER_EVENT
436          
437           if(eventType!=7 && fDebug > 1 )printf("AliCaloTrackReader::FillInputEvent() - DO LED, Event Type <%d>, 8 Calibration \n",  eventType);
438           if(eventType!=8)return kFALSE;
439   }
440                 
441   //In case of analysis of events with jets, skip those with jet pt > 5 pt hard 
442   if(fComparePtHardAndJetPt && GetStack()) {
443     if(!ComparePtHardAndJetPt()) return kFALSE ;
444   }
445
446   //In case of mixing events with other AOD file        
447  // if(fDataType == kAOD && fSecondInputAODTree){
448 //       
449 //    if(fDebug > 1) 
450 //      printf("AliCaloTrackReader::FillInputEvent() - Get event %d from second input AOD file \n", iEntry+fSecondInputFirstEvent);
451 //    if(fSecondInputAODTree->GetEntriesFast() <= iEntry+fSecondInputFirstEvent) {
452 //      if(fSecondInputAODTree->GetEntriesFast() == iEntry+fSecondInputFirstEvent) 
453 //                       printf("AliCaloTrackReader::FillInputEvent() - Skip events from event %d, no more events in second AOD file \n", iEntry);
454 //      return kFALSE;
455 //    }
456 //    
457 //    //Get the Event
458 //    Int_t nbytes = fSecondInputAODTree->GetEvent(iEntry+fSecondInputFirstEvent);
459 //    if ( nbytes == 0 ) {//If nothing in AOD
460 //      printf("AliCaloTrackReader::FillInputEvent() - Nothing in Second AOD input, STOP\n");
461 //      abort() ; 
462 //    }
463 //    
464 //  }
465         
466   //Fill Vertex array
467   FillVertexArray();
468   //Reject events with Z vertex too large, only for SE analysis, if not, cut on the analysis code
469   if(!GetMixedEvent() && TMath::Abs(fVertex[0][2]) > fZvtxCut) return kFALSE;  
470   
471   //------------------------------------------------------
472   //Event rejection depending on vertex, pileup, v0and
473   //------------------------------------------------------
474   if(fDoEventSelection){
475     if(!fCaloFilterPatch){
476       //Do not analyze events with pileup
477       Bool_t bPileup = fInputEvent->IsPileupFromSPD(3, 0.8, 3., 2., 5.); //Default values, if not it does not compile
478       //Bool_t bPileup = event->IsPileupFromSPD(); 
479       if(bPileup) return kFALSE;
480       
481       if(fDoV0ANDEventSelection){
482         Bool_t bV0AND = kTRUE; 
483         AliESDEvent* esd = dynamic_cast<AliESDEvent*> (fInputEvent);
484         if(esd) 
485           bV0AND = fTriggerAnalysis->IsOfflineTriggerFired(esd, AliTriggerAnalysis::kV0AND);
486         //else bV0AND = //FIXME FOR AODs
487         if(!bV0AND) return kFALSE;
488       }
489       
490       if(fUseEventsWithPrimaryVertex && !CheckForPrimaryVertex()) return kFALSE;
491       
492     }//CaloFilter patch
493     else{ 
494       if(fInputEvent->GetNumberOfCaloClusters() > 0) {
495         AliVCluster * calo = fInputEvent->GetCaloCluster(0);
496         if(calo->GetNLabels() == 4){
497           Int_t * selection = calo->GetLabels();
498           Bool_t bPileup = selection[0];
499           if(bPileup) return kFALSE;
500           
501           Bool_t bGoodV = selection[1]; 
502           if(fUseEventsWithPrimaryVertex && !bGoodV) return kFALSE;
503           
504           if(fDoV0ANDEventSelection){
505             Bool_t bV0AND = selection[2]; 
506             if(!bV0AND) return kFALSE;
507           }
508           
509           fTrackMult = selection[3];
510           if(fTrackMult == 0) return kFALSE;
511         } else {
512           //First filtered AODs, track multiplicity stored there.  
513           fTrackMult = (Int_t) ((AliAODHeader*)fInputEvent->GetHeader())->GetCentrality();
514           if(fTrackMult == 0) return kFALSE;          
515         }
516       }//at least one cluster
517       else {
518         //printf("AliCaloTrackReader::FillInputEvent() - No clusters in event\n");
519         //Remove events with  vertex (0,0,0), bad vertex reconstruction
520         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;
521         
522         //First filtered AODs, track multiplicity stored there.  
523         fTrackMult = (Int_t) ((AliAODHeader*)fInputEvent->GetHeader())->GetCentrality();
524         if(fTrackMult == 0) return kFALSE;
525       }// no cluster
526     }// CaloFileter patch
527   }// Event selection
528   //------------------------------------------------------
529
530   //Check if there is a centrality value, PbPb analysis, and if a centrality bin selection is requested
531   //If we need a centrality bin, we select only those events in the corresponding bin.
532   if(GetCentrality() && fCentralityBin[0]>=0 && fCentralityBin[1]>=0 && fCentralityOpt==100){
533     Int_t cen = GetEventCentrality();
534     if(cen > fCentralityBin[1] || cen < fCentralityBin[0]) return kFALSE; //reject events out of bin.
535   }
536   
537   //Fill the arrays with cluster/tracks/cells data
538    if(fFillEMCALCells) 
539     FillInputEMCALCells();
540   if(fFillPHOSCells)  
541     FillInputPHOSCells();
542         
543   if(fFillCTS){   
544     FillInputCTS();
545     //Accept events with at least one track
546     if(fTrackMult == 0 && fDoEventSelection) return kFALSE;
547   }
548   
549   if(fFillEMCAL) 
550     FillInputEMCAL();
551   if(fFillPHOS)  
552     FillInputPHOS();
553
554   FillInputVZERO();
555         
556   return kTRUE ;
557 }
558
559 //__________________________________________________
560 void AliCaloTrackReader::ResetLists() {
561   //  Reset lists, called by the analysis maker 
562
563   if(fCTSTracks)       fCTSTracks     -> Clear();
564   if(fEMCALClusters)   fEMCALClusters   -> Clear("C");
565   if(fPHOSClusters)    fPHOSClusters    -> Clear("C");
566 //  if(fEMCALCells) fEMCALCells -> Clear("");
567 //  if(fPHOSCells)  fPHOSCells  -> Clear("");
568
569   fV0ADC[0] = 0;   fV0ADC[1] = 0; 
570   fV0Mul[0] = 0;   fV0Mul[1] = 0; 
571
572 }
573
574 //____________________________________________________________________________
575 void AliCaloTrackReader::SetInputEvent(AliVEvent* const input)  
576 {
577   fInputEvent  = input;
578   fMixedEvent = dynamic_cast<AliMixedEvent*>(GetInputEvent()) ; 
579   if (fMixedEvent) {
580     fNMixedEvent = fMixedEvent->GetNumberOfEvents() ; 
581   }
582
583   //Delete previous vertex
584   if(fVertex){
585     for (Int_t i = 0; i < fNMixedEvent; i++) {
586       delete [] fVertex[i] ; 
587     }
588     delete [] fVertex ;
589   }
590   
591   fVertex = new Double_t*[fNMixedEvent] ; 
592   for (Int_t i = 0; i < fNMixedEvent; i++) {
593     fVertex[i] = new Double_t[3] ; 
594     fVertex[i][0] = 0.0 ; 
595     fVertex[i][1] = 0.0 ; 
596     fVertex[i][2] = 0.0 ; 
597   }
598 }
599
600 //__________________________________________________
601 Int_t AliCaloTrackReader::GetEventCentrality() const {
602   //Return current event centrality
603   
604   if(GetCentrality()){
605     if(fCentralityOpt==100)     return (Int_t) GetCentrality()->GetCentralityPercentile(fCentralityClass); // 100 bins max
606     else if(fCentralityOpt==10) return GetCentrality()->GetCentralityClass10(fCentralityClass);// 10 bins max
607     else if(fCentralityOpt==20) return GetCentrality()->GetCentralityClass5(fCentralityClass); // 20 bins max
608     else {
609       printf("AliAnaPartCorrBaseClass::Unknown centrality option %d, use 10, 20 or 100\n",fCentralityOpt);
610       return 0;
611     } 
612   }
613   else return 0;
614   
615 }
616
617 //____________________________________________________________________________
618 void AliCaloTrackReader::GetVertex(Double_t vertex[3]) const {
619   //Return vertex position to be used for single event analysis
620   vertex[0]=fVertex[0][0];  
621   vertex[1]=fVertex[0][1];  
622   vertex[2]=fVertex[0][2];
623 }
624
625 //____________________________________________________________________________
626 void AliCaloTrackReader::GetVertex(Double_t vertex[3], const Int_t evtIndex) const {
627   //Return vertex position for mixed event, recover the vertex in a particular event.
628   
629   //Int_t evtIndex = 0; // for single events only one vertex stored in position 0, default value
630   //if (fMixedEvent && clusterID >=0) {
631   //  evtIndex=GetMixedEvent()->EventIndexForCaloCluster(clusterID) ; 
632   //}
633   
634   vertex[0]=fVertex[evtIndex][0];  vertex[1]=fVertex[evtIndex][1];  vertex[2]=fVertex[evtIndex][2];
635   
636 }
637 //
638
639
640 //____________________________________________________________________________
641 void AliCaloTrackReader::FillVertexArray() {
642   
643   //Fill data member with vertex
644   //In case of Mixed event, multiple vertices
645   
646   //Delete previous vertex
647   if(fVertex){
648     for (Int_t i = 0; i < fNMixedEvent; i++) {
649       delete [] fVertex[i] ; 
650     }
651     delete [] fVertex ;  
652   }
653   
654   fVertex = new Double_t*[fNMixedEvent] ; 
655   for (Int_t i = 0; i < fNMixedEvent; i++) {
656     fVertex[i] = new Double_t[3] ; 
657     fVertex[i][0] = 0.0 ; 
658     fVertex[i][1] = 0.0 ; 
659     fVertex[i][2] = 0.0 ; 
660   }          
661   
662   if (!fMixedEvent) { //Single event analysis
663     if(fDataType!=kMC){
664
665       if(fInputEvent->GetPrimaryVertex()){
666         fInputEvent->GetPrimaryVertex()->GetXYZ(fVertex[0]); 
667       }
668       else {
669         printf("AliCaloTrackReader::FillVertexArray() - NULL primary vertex\n");
670         fVertex[0][0]=0.;   fVertex[0][1]=0.;   fVertex[0][2]=0.;
671       }//Primary vertex pointer do not exist
672       
673     } else {//MC read event 
674       fVertex[0][0]=0.;   fVertex[0][1]=0.;   fVertex[0][2]=0.;
675     }
676       
677     if(fDebug > 1)
678       printf("AliCaloTrackReader::FillVertexArray() - Single Event Vertex : %f,%f,%f\n",fVertex[0][0],fVertex[0][1],fVertex[0][2]);
679
680   } else { // MultiEvent analysis
681     for (Int_t iev = 0; iev < fNMixedEvent; iev++) {
682       if (fMixedEvent->GetVertexOfEvent(iev))
683         fMixedEvent->GetVertexOfEvent(iev)->GetXYZ(fVertex[iev]);
684       else { // no vertex found !!!!
685         AliWarning("No vertex found");
686       }
687
688       if(fDebug > 1)
689         printf("AliCaloTrackReader::FillVertexArray() - Multi Event %d Vertex : %f,%f,%f\n",iev,fVertex[iev][0],fVertex[iev][1],fVertex[iev][2]);
690
691     }
692   }
693   
694 }
695
696 //____________________________________________________________________________
697 void AliCaloTrackReader::FillInputCTS() {
698   //Return array with Central Tracking System (CTS) tracks
699   
700   if(fDebug > 2 ) printf("AliCaloTrackReader::FillInputCTS()\n");
701   
702   Int_t nTracks   = fInputEvent->GetNumberOfTracks() ;
703   Double_t p[3];
704   fTrackMult = 0;
705   Int_t nstatus = 0;
706   for (Int_t itrack =  0; itrack <  nTracks; itrack++) {////////////// track loop
707     AliVTrack * track = (AliVTrack*)fInputEvent->GetTrack(itrack) ; // retrieve track from esd
708
709     //Select tracks under certain conditions, TPCrefit, ITSrefit ... check the set bits
710     if (fTrackStatus && !((track->GetStatus() & fTrackStatus) == fTrackStatus)) 
711       continue ;
712     
713     nstatus++;
714     
715     if(fDataType==kESD && !fESDtrackCuts->AcceptTrack((AliESDtrack*)track)) continue;
716     
717     // Track filter selection
718     //if (fTrackFilter) {
719           //  selectInfo = fTrackFilter->IsSelected(esdTrack);
720           //  if (!selectInfo && !(esd->GetPrimaryVertex())->UsesTrack(esdTrack->GetID())) continue;
721    // }
722     
723     //Count the tracks in eta < 0.9
724     //printf("Eta %f cut  %f\n",TMath::Abs(track->Eta()),fTrackMultEtaCut);
725     if(TMath::Abs(track->Eta())< fTrackMultEtaCut) fTrackMult++;
726     
727     track->GetPxPyPz(p) ;
728     TLorentzVector momentum(p[0],p[1],p[2],0);
729     
730     if(fCTSPtMin < momentum.Pt() && fCTSPtMax > momentum.Pt()){
731       
732       if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"CTS")) 
733         continue;
734       
735       if(fDebug > 2 && momentum.Pt() > 0.1) 
736         printf("AliCaloTrackReader::FillInputCTS() - Selected tracks E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
737                momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
738       
739       if (fMixedEvent) {
740         track->SetID(itrack);
741       }
742       
743       fCTSTracks->Add(track);        
744        
745     }//Pt and Fiducial cut passed. 
746   }// track loop
747         
748   //fCTSTracksNormalInputEntries = fCTSTracks->GetEntriesFast();
749   if(fDebug > 1) 
750     printf("AliCaloTrackReader::FillInputCTS()   - aod entries %d, input tracks %d, pass status %d, multipliticy %d\n", fCTSTracks->GetEntriesFast(), nTracks, nstatus, fTrackMult);//fCTSTracksNormalInputEntries);
751   
752     //  //If second input event available, add the clusters.
753     //  if(fSecondInputAODTree && fSecondInputAODEvent){
754     //    nTracks   = fSecondInputAODEvent->GetNumberOfTracks() ;
755     //    if(fDebug > 1) printf("AliCaloTrackReader::FillInputCTS()   - Add second input tracks, entries %d\n", nTracks) ;
756     //    for (Int_t itrack =  0; itrack <  nTracks; itrack++) {////////////// track loop
757     //            AliAODTrack * track = ((AliAODEvent*)fSecondInputAODEvent)->GetTrack(itrack) ; // retrieve track from esd
758     //            
759     //            //Select tracks under certain conditions, TPCrefit, ITSrefit ... check the set bits
760     //            if (fTrackStatus && !((track->GetStatus() & fTrackStatus) == fTrackStatus)) continue ;
761     //            
762     //            track->GetPxPyPz(p) ;
763     //            TLorentzVector momentum(p[0],p[1],p[2],0);
764     //            
765     //            if(fCTSPtMin < momentum.Pt()){
766     //
767     //                    if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"CTS")) continue;
768     //
769     //                    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",
770     //                                                                 momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
771     //                    
772     //                    fCTSTracks->Add(track);
773     //                    
774     //            }//Pt and Fiducial cut passed. 
775     //    }// track loop
776     //    
777     //    if(fDebug > 1) printf("AliCaloTrackReader::FillInputCTS()   - aod normal entries %d, after second input %d\n", fCTSTracksNormalInputEntries, fCTSTracks->GetEntriesFast());
778     //  }       //second input loop
779     //  
780 }
781
782 //____________________________________________________________________________
783 void AliCaloTrackReader::FillInputEMCALAlgorithm(AliVCluster * clus, const Int_t iclus) {
784   //Fill the EMCAL data in the array, do it
785   
786   Int_t vindex = 0 ;  
787   if (fMixedEvent) 
788     vindex = fMixedEvent->EventIndexForCaloCluster(iclus);
789   
790   
791   //Reject clusters with bad channels, close to borders and exotic;
792   if(!GetCaloUtils()->GetEMCALRecoUtils()->IsGoodCluster(clus,GetCaloUtils()->GetEMCALGeometry(),GetEMCALCells())) return;
793 //  //Check if the cluster contains any bad channel and if close to calorimeter borders
794 //  if(GetCaloUtils()->ClusterContainsBadChannel("EMCAL",clus->GetCellsAbsId(), clus->GetNCells())) 
795 //    return;
796 //  if(!GetCaloUtils()->CheckCellFiducialRegion(clus, (AliVCaloCells*)fInputEvent->GetEMCALCells(), fInputEvent, vindex)) 
797 //    return;
798 //  
799   
800   if(fSelectEmbeddedClusters){
801     if(clus->GetNLabels()==0 || clus->GetLabel() < 0) return;
802     //else printf("Embedded cluster,  %d, n label %d label %d  \n",iclus,clus->GetNLabels(),clus->GetLabel());
803   }
804   
805   TLorentzVector momentum ;
806   
807   clus->GetMomentum(momentum, fVertex[vindex]);      
808   
809   if(fEMCALPtMin < momentum.E() && fEMCALPtMax > momentum.E()){
810     
811     if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"EMCAL")) 
812       return;
813     
814     if(fDebug > 2 && momentum.E() > 0.1) 
815       printf("AliCaloTrackReader::FillInputEMCAL() - Selected clusters E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
816              momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
817     
818     //Float_t pos[3];
819     //clus->GetPosition(pos);
820     //printf("Before Corrections: e %f, x %f, y %f, z %f\n",clus->E(),pos[0],pos[1],pos[2]);
821     
822     //Recalibrate the cluster energy 
823     if(GetCaloUtils()->IsRecalibrationOn()) {
824       Float_t energy = GetCaloUtils()->RecalibrateClusterEnergy(clus, GetEMCALCells());
825       clus->SetE(energy);
826       //printf("Recalibrated Energy %f\n",clus->E());  
827       GetCaloUtils()->RecalculateClusterShowerShapeParameters(GetEMCALCells(),clus);
828       GetCaloUtils()->RecalculateClusterPID(clus);
829       
830     }
831     
832     //Recalculate distance to bad channels, if new list of bad channels provided
833     GetCaloUtils()->RecalculateClusterDistanceToBadChannel(GetEMCALCells(),clus);
834     
835     //Recalculate cluster position
836     if(GetCaloUtils()->IsRecalculationOfClusterPositionOn()){
837       GetCaloUtils()->RecalculateClusterPosition(GetEMCALCells(),clus); 
838       //clus->GetPosition(pos);
839       //printf("After  Corrections: e %f, x %f, y %f, z %f\n",clus->E(),pos[0],pos[1],pos[2]);
840     }
841     
842     //Correct non linearity
843     if(GetCaloUtils()->IsCorrectionOfClusterEnergyOn()){
844       GetCaloUtils()->CorrectClusterEnergy(clus) ;
845       //printf("Linearity Corrected Energy %f\n",clus->E());  
846     }          
847     
848     //In case of MC analysis, to match resolution/calibration in real data
849     if(fSmearClusterEnergy){
850       Float_t energy    = clus->E();
851       Float_t rdmEnergy = fRandom.Gaus(energy,fSmearClusterParam[0]*TMath::Sqrt(energy)+
852                                        fSmearClusterParam[1]*energy+fSmearClusterParam[2]);
853       clus->SetE(rdmEnergy);
854       if(fDebug > 2) printf("\t Energy %f, smeared %f\n", energy, clus->E());
855     }
856     
857     if (fMixedEvent) 
858       clus->SetID(iclus) ; 
859     
860     fEMCALClusters->Add(clus);  
861   }
862 }
863
864 //____________________________________________________________________________
865 void AliCaloTrackReader::FillInputEMCAL() {
866   //Return array with EMCAL clusters in aod format
867   
868   if(fDebug > 2 ) printf("AliCaloTrackReader::FillInputEMCAL()\n");
869   
870   //Loop to select clusters in fiducial cut and fill container with aodClusters
871   if(fEMCALClustersListName==""){
872     Int_t nclusters = fInputEvent->GetNumberOfCaloClusters();
873     for (Int_t iclus =  0; iclus <  nclusters; iclus++) {
874       AliVCluster * clus = 0;
875       if ( (clus = fInputEvent->GetCaloCluster(iclus)) ) {
876         if (IsEMCALCluster(clus)){          
877           FillInputEMCALAlgorithm(clus, iclus);
878         }//EMCAL cluster
879       }// cluster exists
880     }// cluster loop
881     
882     //Recalculate track matching
883     GetCaloUtils()->RecalculateClusterTrackMatching(fInputEvent);
884     
885   }//Get the clusters from the input event
886   else {
887     TClonesArray * clusterList = 0x0; 
888     if(fOutputEvent) clusterList = dynamic_cast<TClonesArray*> (fOutputEvent->FindListObject(fEMCALClustersListName));
889     if(!clusterList){
890       //printf("AliCaloTrackReader::FillInputEMCAL() - Wrong name of list with clusters? Try input event <%s>\n",fEMCALClustersListName.Data());
891       //List not in output event, try input event
892       clusterList = dynamic_cast<TClonesArray*> (fInputEvent->FindListObject(fEMCALClustersListName));
893       if(!clusterList){
894         printf("AliCaloTrackReader::FillInputEMCAL() - Wrong name of list with clusters?  <%s>\n",fEMCALClustersListName.Data());
895         return;
896       }
897     }
898     Int_t nclusters = clusterList->GetEntriesFast();
899     for (Int_t iclus =  0; iclus <  nclusters; iclus++) {
900       AliVCluster * clus = dynamic_cast<AliVCluster*> (clusterList->At(iclus));
901       //printf("E %f\n",clus->E());
902       if (clus) FillInputEMCALAlgorithm(clus, iclus);
903       else printf("AliCaloTrackReader::FillInputEMCAL() - Null cluster in list!\n");
904       
905     }// cluster loop
906     
907     //Recalculate track matching, not necessary, already done in the reclusterization task
908     //GetCaloUtils()->RecalculateClusterTrackMatching(fInputEvent,clusterList);
909     
910   }
911     
912   //fEMCALClustersNormalInputEntries = fEMCALClusters->GetEntriesFast();
913   if(fDebug > 1) printf("AliCaloTrackReader::FillInputEMCAL() - aod entries %d\n",  fEMCALClusters->GetEntriesFast());//fEMCALClustersNormalInputEntries);
914   
915     //If second input event available, add the clusters.
916     //  if(fSecondInputAODTree && fSecondInputAODEvent){
917     //    GetSecondInputAODVertex(v);
918     //    nclusters = ((AliAODEvent*)fSecondInputAODEvent)->GetNumberOfCaloClusters();
919     //    if(fDebug > 1) printf("AliCaloTrackReader::FillInputEMCAL() - Add second input clusters, entries %d\n", nclusters) ;
920     //          for (Int_t iclus =  0; iclus < nclusters; iclus++) {
921     //                  AliAODCaloCluster * clus = 0;
922     //                  if ( (clus = ((AliAODEvent*)fSecondInputAODEvent)->GetCaloCluster(iclus)) ) {
923     //                          if (clus->IsEMCAL()){
924     //                                  TLorentzVector momentum ;
925     //                                  clus->GetMomentum(momentum, v);      
926     //                                  
927     //                                  if(fEMCALPtMin < momentum.Pt()){
928     //
929     //                                          if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"EMCAL")) continue;
930     //
931     //                                          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",
932     //                                                                                                                                  momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
933     //                                    fEMCALClusters->Add(clus);    
934     //                                  }//Pt and Fiducial cut passed.
935     //                          }//EMCAL cluster
936     //                  }// cluster exists
937     //          }// cluster loop
938     //          
939     //    if(fDebug > 1) printf("AliCaloTrackReader::FillInputEMCAL() - aod normal entries %d, after second input %d\n", fEMCALClustersNormalInputEntries, fEMCALClusters->GetEntriesFast());
940     //
941     //  } //second input loop
942 }
943
944   //____________________________________________________________________________
945 void AliCaloTrackReader::FillInputPHOS() {
946   //Return array with PHOS clusters in aod format
947   
948   if(fDebug > 2 ) printf("AliCaloTrackReader::FillInputPHOS()\n");
949           
950   //Loop to select clusters in fiducial cut and fill container with aodClusters
951   Int_t nclusters = fInputEvent->GetNumberOfCaloClusters();
952   for (Int_t iclus = 0; iclus < nclusters; iclus++) {
953     AliVCluster * clus = 0;
954     if ( (clus = fInputEvent->GetCaloCluster(iclus)) ) {
955       if (IsPHOSCluster(clus)){
956         //Check if the cluster contains any bad channel and if close to calorimeter borders
957         Int_t vindex = 0 ;  
958         if (fMixedEvent) 
959           vindex = fMixedEvent->EventIndexForCaloCluster(iclus);
960         if( GetCaloUtils()->ClusterContainsBadChannel("PHOS",clus->GetCellsAbsId(), clus->GetNCells())) 
961           continue;
962         if(!GetCaloUtils()->CheckCellFiducialRegion(clus, fInputEvent->GetPHOSCells(), fInputEvent, vindex)) 
963           continue;
964         
965         TLorentzVector momentum ;
966         
967         clus->GetMomentum(momentum, fVertex[vindex]);      
968         
969         if(fPHOSPtMin < momentum.E() && fPHOSPtMax > momentum.E()){
970           
971           if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"PHOS")) 
972             continue;
973           
974           if(fDebug > 2 && momentum.E() > 0.1) 
975             printf("AliCaloTrackReader::FillInputPHOS() - Selected clusters E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
976                    momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
977           
978             //Recalibrate the cluster energy 
979           if(GetCaloUtils()->IsRecalibrationOn()) {
980             Float_t energy = GetCaloUtils()->RecalibrateClusterEnergy(clus, (AliAODCaloCells*)GetPHOSCells());
981             clus->SetE(energy);
982           }
983           
984           if (fMixedEvent) {
985             clus->SetID(iclus) ; 
986           }              
987           
988           fPHOSClusters->Add(clus);     
989           
990         }//Pt and Fiducial cut passed.
991       }//PHOS cluster
992     }//cluster exists
993   }//esd cluster loop
994   
995   //fPHOSClustersNormalInputEntries = fPHOSClusters->GetEntriesFast() ;
996   if(fDebug > 1) printf("AliCaloTrackReader::FillInputPHOS()  - aod entries %d\n",  fPHOSClusters->GetEntriesFast());//fPHOSClustersNormalInputEntries);
997   
998     //If second input event available, add the clusters.
999     //  if(fSecondInputAODTree && fSecondInputAODEvent){  
1000     //    GetSecondInputAODVertex(v);
1001     //    nclusters = ((AliAODEvent*)fSecondInputAODEvent)->GetNumberOfCaloClusters();
1002     //    if(fDebug > 1) printf("AliCaloTrackReader::FillInputPHOS()  - Add second input clusters, entries %d\n", nclusters);
1003     //          for (Int_t iclus =  0; iclus < nclusters; iclus++) {
1004     //                  AliAODCaloCluster * clus = 0;
1005     //                  if ( (clus = ((AliAODEvent*)fSecondInputAODEvent)->GetCaloCluster(iclus)) ) {
1006     //                          if (clus->IsPHOS()){
1007     //                                  TLorentzVector momentum ;
1008     //                                  clus->GetMomentum(momentum, v);      
1009     //                                  
1010     //                                  if(fPHOSPtMin < momentum.Pt()){
1011     //
1012     //                                          if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"PHOS")) continue;
1013     //
1014     //                                          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",
1015     //                                                                                                                                  momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
1016     //                                          fPHOSClusters->Add(clus);       
1017     //                                  }//Pt and Fiducial cut passed.
1018     //                          }//PHOS cluster
1019     //                  }// cluster exists
1020     //          }// cluster loop
1021     //          if(fDebug > 1) printf("AliCaloTrackReader::FillInputPHOS()  - aod normal entries %d, after second input %d\n", fPHOSClustersNormalInputEntries, fPHOSClusters->GetEntriesFast());
1022     //  }       //second input loop
1023   
1024 }
1025
1026 //____________________________________________________________________________
1027 void AliCaloTrackReader::FillInputEMCALCells() {
1028     //Return array with EMCAL cells in aod format
1029   
1030   fEMCALCells = fInputEvent->GetEMCALCells(); 
1031   
1032 }
1033
1034 //____________________________________________________________________________
1035 void AliCaloTrackReader::FillInputPHOSCells() {
1036     //Return array with PHOS cells in aod format
1037   
1038   fPHOSCells = fInputEvent->GetPHOSCells(); 
1039   
1040 }
1041
1042 //____________________________________________________________________________
1043 void AliCaloTrackReader::FillInputVZERO(){
1044   //Fill VZERO information in data member, add all the channels information.
1045   AliVVZERO* v0 = fInputEvent->GetVZEROData();
1046   //printf("Init V0: ADC (%d,%d), Multiplicity (%d,%d) \n",fV0ADC[0],fV0ADC[1],fV0Mul[0],fV0Mul[1]);
1047   
1048   if (v0) 
1049   {
1050     AliESDVZERO* esdV0 = dynamic_cast<AliESDVZERO*> (v0);
1051     for (Int_t i = 0; i < 32; i++)
1052     {
1053       if(esdV0){//Only available in ESDs
1054         fV0ADC[0] += (Int_t)esdV0->GetAdcV0C(i);
1055         fV0ADC[1] += (Int_t)esdV0->GetAdcV0A(i);
1056       }
1057       fV0Mul[0] += (Int_t)v0->GetMultiplicityV0C(i);
1058       fV0Mul[1] += (Int_t)v0->GetMultiplicityV0A(i);
1059     }
1060     if(fDebug > 0)
1061       printf("V0: ADC (%d,%d), Multiplicity (%d,%d) \n",fV0ADC[0],fV0ADC[1],fV0Mul[0],fV0Mul[1]);
1062   }
1063   else
1064   {
1065     if(fDebug > 0)
1066       printf("Cannot retrieve V0 ESD! Run w/ null V0 charges\n ");
1067   }
1068 }
1069
1070
1071 //____________________________________________________________________________
1072 Bool_t AliCaloTrackReader::IsEMCALCluster(AliVCluster* cluster) const {
1073   // Check if it is a cluster from EMCAL. For old AODs cluster type has
1074   // different number and need to patch here
1075     
1076   if(fDataType==kAOD && fOldAOD)
1077   {
1078     if (cluster->GetType() == 2) return kTRUE;
1079     else                         return kFALSE;
1080   }
1081   else 
1082   {
1083     return cluster->IsEMCAL();
1084   }
1085
1086 }
1087
1088 //____________________________________________________________________________
1089 Bool_t AliCaloTrackReader::IsPHOSCluster(AliVCluster * cluster) const {
1090   //Check if it is a cluster from PHOS.For old AODs cluster type has
1091   // different number and need to patch here
1092   
1093   if(fDataType==kAOD && fOldAOD)
1094   {
1095     Int_t type = cluster->GetType();
1096     if (type == 0 || type == 1) return kTRUE;
1097     else                        return kFALSE;
1098   }
1099   else 
1100   {
1101     return cluster->IsPHOS();
1102   }
1103   
1104 }
1105
1106 //____________________________________________________________________________
1107 Bool_t AliCaloTrackReader::CheckForPrimaryVertex(){
1108   //Check if the vertex was well reconstructed, copy from V0Reader of conversion group
1109   //Only for ESDs ...
1110   
1111   AliESDEvent * event = dynamic_cast<AliESDEvent*> (fInputEvent);
1112   if(!event) return kTRUE;
1113   
1114   if(event->GetPrimaryVertexTracks()->GetNContributors() > 0) {
1115     return kTRUE;
1116   }
1117   
1118   if(event->GetPrimaryVertexTracks()->GetNContributors() < 1) {
1119     // SPD vertex
1120     if(event->GetPrimaryVertexSPD()->GetNContributors() > 0) {
1121       //cout<<"spd vertex type::"<< fESDEvent->GetPrimaryVertex()->GetName() << endl;
1122       return kTRUE;
1123       
1124     }
1125     if(event->GetPrimaryVertexSPD()->GetNContributors() < 1) {
1126       //      cout<<"bad vertex type::"<< event->GetPrimaryVertex()->GetName() << endl;
1127       return kFALSE;
1128     }
1129   }
1130   
1131   return kFALSE;  
1132   
1133 }
1134
1135
1136