]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/PartCorrBase/AliCaloTrackReader.cxx
get the v0 signal and multiplicity, accessible for different analysis, add histograms...
[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
35 //---- ANALYSIS system ----
36 #include "AliCaloTrackReader.h"
37 #include "AliMCEvent.h"
38 #include "AliAODMCHeader.h"
39 #include "AliGenPythiaEventHeader.h"
40 #include "AliVEvent.h"
41 #include "AliAODEvent.h"
42 #include "AliVTrack.h"
43 #include "AliVParticle.h"
44 #include "AliMixedEvent.h"
45 #include "AliESDtrack.h"
46 #include "AliEMCALRecoUtils.h"
47 #include "AliESDtrackCuts.h"
48
49 ClassImp(AliCaloTrackReader)
50   
51   
52 //____________________________________________________________________________
53   AliCaloTrackReader::AliCaloTrackReader() : 
54     TObject(), fEventNumber(-1), fCurrentFileName(""),fDataType(0), fDebug(0), 
55     fFiducialCut(0x0), fCheckFidCut(kFALSE), fComparePtHardAndJetPt(kFALSE), fPtHardAndJetPtFactor(7),
56     fCTSPtMin(0), fEMCALPtMin(0),fPHOSPtMin(0), fAODBranchList(new TList ),
57     fAODCTS(new TObjArray()), fAODEMCAL(new TObjArray()), fAODPHOS(new TObjArray()),
58     fEMCALCells(0x0), fPHOSCells(0x0),
59     fInputEvent(0x0), fOutputEvent(0x0),fMC(0x0),
60     fFillCTS(0),fFillEMCAL(0),fFillPHOS(0),
61     fFillEMCALCells(0),fFillPHOSCells(0), 
62 //    fSecondInputAODTree(0x0), fSecondInputAODEvent(0x0),
63 //    fSecondInputFileName(""),fSecondInputFirstEvent(0), 
64 //    fAODCTSNormalInputEntries(0), fAODEMCALNormalInputEntries(0), 
65 //    fAODPHOSNormalInputEntries(0), 
66     fTrackStatus(0),   fESDtrackCuts(0), fTrackMult(0), fTrackMultEtaCut(0.9),
67     fReadStack(kFALSE), fReadAODMCParticles(kFALSE), 
68     fDeltaAODFileName("deltaAODPartCorr.root"),fFiredTriggerClassName(""),
69     fAnaLED(kFALSE),fTaskName(""),fCaloUtils(0x0), 
70     fMixedEvent(NULL), fNMixedEvent(1), fVertex(NULL), 
71     fWriteOutputDeltaAOD(kFALSE),fOldAOD(kFALSE){
72   //Ctor
73   
74   //Initialize parameters
75   InitParameters();
76 }
77
78 //_________________________________
79 AliCaloTrackReader::~AliCaloTrackReader() {
80   //Dtor
81   
82   if(fFiducialCut) delete fFiducialCut ;
83         
84   if(fAODBranchList){
85     fAODBranchList->Delete();
86     delete fAODBranchList ;
87   }  
88   
89   if(fAODCTS){
90     if(fDataType!=kMC)fAODCTS->Clear() ; 
91     else              fAODCTS->Delete() ; 
92     delete fAODCTS ;
93   }
94   
95   if(fAODEMCAL){
96     if(fDataType!=kMC)fAODEMCAL->Clear() ; 
97     else              fAODEMCAL->Delete() ; 
98     delete fAODEMCAL ;
99   }
100   
101   if(fAODPHOS){
102     if(fDataType!=kMC)fAODPHOS->Clear() ; 
103     else              fAODPHOS->Delete() ; 
104     delete fAODPHOS ;
105   }
106   
107   if(fEMCALCells){
108     fEMCALCells->Clear() ; 
109     delete fEMCALCells ;
110   }
111   
112   if(fPHOSCells){
113     fPHOSCells->Clear() ; 
114     delete fPHOSCells ;
115   }
116
117   if(fVertex){
118     for (Int_t i = 0; i < fNMixedEvent; i++) {
119       delete [] fVertex[i] ;
120
121     }
122     delete [] fVertex ;
123         }
124
125   if(fESDtrackCuts)   delete fESDtrackCuts;
126   
127 //  Pointers not owned, done by the analysis frame
128 //  if(fInputEvent)  delete fInputEvent ;
129 //  if(fOutputEvent) delete fOutputEvent ;
130 //  if(fMC)          delete fMC ;  
131         
132 //  if(fSecondInputAODTree){
133 //    fSecondInputAODTree->Clear();
134 //    delete fSecondInputAODTree;
135 //  }
136 //      
137 //  if(fSecondInputAODEvent) delete fSecondInputAODEvent ;
138         
139   //  Pointer not owned, deleted by maker
140   //if (fCaloUtils) delete fCaloUtils ;
141
142 }
143
144 //_________________________________________________________________________
145 Bool_t AliCaloTrackReader::ComparePtHardAndJetPt(){
146         // Check the event, if the requested ptHard is much larger than the jet pT, then there is a problem.
147         // Only for PYTHIA.
148         if(!fReadStack) return kTRUE; //Information not filtered to AOD
149         
150         if(!strcmp(GetGenEventHeader()->ClassName(), "AliGenPythiaEventHeader")){
151                 TParticle * jet =  0;
152                 AliGenPythiaEventHeader* pygeh= (AliGenPythiaEventHeader*) GetGenEventHeader();
153                 Int_t nTriggerJets =  pygeh->NTriggerJets();
154                 Float_t ptHard = pygeh->GetPtHard();
155                 
156                 //if(fDebug > 1) printf("AliMCAnalysisUtils::PythiaEventHeader: Njets: %d, pT Hard %f\n",nTriggerJets, ptHard);
157     Float_t tmpjet[]={0,0,0,0};
158                 for(Int_t ijet = 0; ijet< nTriggerJets; ijet++){
159                         pygeh->TriggerJet(ijet, tmpjet);
160                         jet = new TParticle(94, 21, -1, -1, -1, -1, tmpjet[0],tmpjet[1],tmpjet[2],tmpjet[3], 0,0,0,0);
161                         //Compare jet pT and pt Hard
162                         //if(fDebug > 1) printf("AliMCAnalysisUtils:: %d pycell jet pT %f\n",ijet, jet->Pt());
163                         if(jet->Pt() > fPtHardAndJetPtFactor * ptHard) {
164                                 printf("AliMCAnalysisUtils::PythiaEventHeader: Njets: %d, pT Hard %2.2f, pycell jet pT %2.2f, rejection factor %1.1f\n",
165                                            nTriggerJets, ptHard, jet->Pt(), fPtHardAndJetPtFactor);
166                                 return kFALSE;
167                         }
168                 }
169     if(jet) delete jet; 
170         }
171           
172         return kTRUE ;
173         
174 }
175
176 //____________________________________________________________________________
177 AliStack* AliCaloTrackReader::GetStack() const {
178   //Return pointer to stack
179   if(fMC)
180     return fMC->Stack();
181   else{
182     if(fDebug > 1) printf("AliCaloTrackReader::GetStack() - Stack is not available\n"); 
183     return 0x0 ;
184   }
185 }
186
187 //____________________________________________________________________________
188 AliHeader* AliCaloTrackReader::GetHeader() const {
189   //Return pointer to header
190   if(fMC)
191     return fMC->Header();
192   else{
193     printf("AliCaloTrackReader::Header is not available\n"); 
194     return 0x0 ;
195   }
196 }
197 //____________________________________________________________________________
198 AliGenEventHeader* AliCaloTrackReader::GetGenEventHeader() const {
199   //Return pointer to Generated event header
200   if(fMC)
201     return fMC->GenEventHeader();
202   else{
203     printf("AliCaloTrackReader::GenEventHeader is not available\n"); 
204     return 0x0 ;
205   }
206 }
207
208 //____________________________________________________________________________
209 TClonesArray* AliCaloTrackReader::GetAODMCParticles(Int_t input) const {
210   //Return list of particles in AOD. Do it for the corresponding input event.
211   
212   TClonesArray * rv = NULL ; 
213   if(fDataType == kAOD){
214  
215     if(input == 0){
216       //Normal input AOD
217       AliAODEvent * evt = dynamic_cast<AliAODEvent*> (fInputEvent) ;
218       if(evt)
219         rv = (TClonesArray*)evt->FindListObject("mcparticles");
220       else  
221         printf("AliCaloTrackReader::GetAODMCParticles() - wrong AOD input index? %d, or non existing tree? \n",input); 
222       
223     } //else if(input == 1 && fSecondInputAODEvent){ //Second input AOD   
224 //      rv = (TClonesArray*) fSecondInputAODEvent->FindListObject("mcparticles");       
225 //    } 
226     
227   } else {
228       printf("AliCaloTrackReader::GetAODMCParticles() - Input are not AODs\n"); 
229   }
230   
231   return rv ; 
232 }
233
234 //____________________________________________________________________________
235 AliAODMCHeader* AliCaloTrackReader::GetAODMCHeader(Int_t input) const {
236         //Return MC header in AOD. Do it for the corresponding input event.
237   AliAODMCHeader *mch = NULL;
238         if(fDataType == kAOD){
239                 //Normal input AOD
240                 if(input == 0) {
241       mch = (AliAODMCHeader*)((AliAODEvent*)fInputEvent)->FindListObject("mcheader");
242     }
243 //              //Second input AOD
244 //              else if(input == 1){ 
245 //       mch = (AliAODMCHeader*) fSecondInputAODEvent->FindListObject("mcheader");
246 //  }
247                 else {
248                         printf("AliCaloTrackReader::GetAODMCHeader() - wrong AOD input index, %d\n",input);
249                 }
250         }
251         else {
252                 printf("AliCaloTrackReader::GetAODMCHeader() - Input are not AODs\n");
253         }
254   
255   return mch;
256 }
257
258 //_______________________________________________________________
259 void AliCaloTrackReader::Init()
260 {
261         //Init reader. Method to be called in AliAnaPartCorrMaker
262         
263         //Get the file with second input events if the filename is given
264         //Get the tree and connect the AODEvent. Only with AODs
265
266         if(fReadStack && fReadAODMCParticles){
267                 printf("AliCaloTrackReader::Init() - Cannot access stack and mcparticles at the same time, change them \n");
268                 fReadStack = kFALSE;
269                 fReadAODMCParticles = kFALSE;
270         }
271         
272 //      if(fSecondInputFileName!=""){
273 //              if(fDataType == kAOD){
274 //                      TFile * input2   = new TFile(fSecondInputFileName,"read");
275 //                      printf("AliCaloTrackReader::Init() - Second input file opened: %s, size %d \n", input2->GetName(), (Int_t) input2->GetSize());
276 //                      fSecondInputAODTree = (TTree*) input2->Get("aodTree");
277 //                      if(fSecondInputAODTree) printf("AliCaloTrackReader::Init() - Second input tree opened: %s, entries %d \n", 
278 //                                                                                 fSecondInputAODTree->GetName(), (Int_t) fSecondInputAODTree->GetEntries());
279 //                      else{
280 //                       printf("AliCaloTrackReader::Init() - Second input tree not available, STOP \n");
281 //                       abort();
282 //                      }
283 //                      fSecondInputAODEvent = new AliAODEvent;
284 //                      fSecondInputAODEvent->ReadFromTree(fSecondInputAODTree);
285 //                      if(fSecondInputFirstEvent >= fSecondInputAODTree->GetEntriesFast()){
286 //                              printf("AliCaloTrackReader::Init() - Requested first event of second input %d, is larger than number of events %d, STOP\n", 
287 //                                         fSecondInputFirstEvent, (Int_t) fSecondInputAODTree->GetEntriesFast());
288 //                              abort();
289 //                      }
290 //              }
291 //              else printf("AliCaloTrackReader::Init() - Second input not added, reader is not AOD\n");
292 //      }
293 }
294
295 //_______________________________________________________________
296 void AliCaloTrackReader::InitParameters()
297 {
298   //Initialize the parameters of the analysis.
299   fDataType   = kESD ;
300   fCTSPtMin   = 0.2 ;
301   fEMCALPtMin = 0.2 ;
302   fPHOSPtMin  = 0.2 ;
303
304   //Do not filter the detectors input by default.
305   fFillEMCAL      = kFALSE;
306   fFillPHOS       = kFALSE;
307   fFillCTS        = kFALSE;
308   fFillEMCALCells = kFALSE;
309   fFillPHOSCells  = kFALSE;
310
311   //fSecondInputFileName   = "" ;
312   //fSecondInputFirstEvent = 0 ;
313   fReadStack             = kFALSE; // Check in the constructor of the other readers if it was set or in the configuration file
314   fReadAODMCParticles    = kFALSE; // Check in the constructor of the other readers if it was set or in the configuration file
315   fDeltaAODFileName      = "deltaAODPartCorr.root";
316   fFiredTriggerClassName      = "";
317                 
318   fAnaLED = kFALSE;
319   
320   //We want tracks fitted in the detectors:
321   //fTrackStatus=AliESDtrack::kTPCrefit;
322   //fTrackStatus|=AliESDtrack::kITSrefit;  
323   
324   fESDtrackCuts = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010();
325
326   fV0ADC[0] = 0;   fV0ADC[1] = 0; 
327   fV0Mul[0] = 0;   fV0Mul[1] = 0; 
328
329 }
330
331 //________________________________________________________________
332 void AliCaloTrackReader::Print(const Option_t * opt) const
333 {
334
335   //Print some relevant parameters set for the analysis
336   if(! opt)
337     return;
338
339   printf("***** Print: %s %s ******\n", GetName(), GetTitle() ) ;
340   printf("Task name      : %s\n", fTaskName.Data()) ;
341   printf("Data type      : %d\n", fDataType) ;
342   printf("CTS Min pT     : %2.1f GeV/c\n", fCTSPtMin) ;
343   printf("EMCAL Min pT   : %2.1f GeV/c\n", fEMCALPtMin) ;
344   printf("PHOS Min pT    : %2.1f GeV/c\n", fPHOSPtMin) ;
345   printf("Use CTS         =     %d\n", fFillCTS) ;
346   printf("Use EMCAL       =     %d\n", fFillEMCAL) ;
347   printf("Use PHOS        =     %d\n", fFillPHOS) ;
348   printf("Use EMCAL Cells =     %d\n", fFillEMCALCells) ;
349   printf("Use PHOS  Cells =     %d\n", fFillPHOSCells) ;
350   printf("Track status    =     %d\n", (Int_t) fTrackStatus) ;
351   printf("Track Mult Eta Cut =  %d\n", (Int_t) fTrackMultEtaCut) ;
352   printf("Write delta AOD =     %d\n", fWriteOutputDeltaAOD) ;
353
354   if(fComparePtHardAndJetPt)
355           printf("Compare jet pt and pt hard to accept event, factor = %2.2f",fPtHardAndJetPtFactor);
356         
357 //  if(fSecondInputFileName!="") {
358 //        printf("Second Input File Name     =     %s\n", fSecondInputFileName.Data()) ;
359 //        printf("Second Input First Event   =     %d\n", fSecondInputFirstEvent) ;
360 //  }
361         
362   printf("Read Kine from, stack? %d, AOD ? %d \n", fReadStack, fReadAODMCParticles) ;
363   printf("Delta AOD File Name =     %s\n", fDeltaAODFileName.Data()) ;
364   printf("    \n") ;
365
366
367 //___________________________________________________
368 Bool_t AliCaloTrackReader::FillInputEvent(const Int_t iEntry, const char * currentFileName) {
369   //Fill the event counter and input lists that are needed, called by the analysis maker.
370
371   fEventNumber = iEntry;
372   fCurrentFileName = TString(currentFileName);
373   if(!fInputEvent) {
374           if(fDebug >= 0) printf("AliCaloTrackReader::FillInputEvent() - Input event not available, skip event analysis\n");
375           return kFALSE;
376   }
377   //Select events only fired by a certain trigger configuration if it is provided
378   Int_t eventType = 0;
379   if(fInputEvent->GetHeader())
380           eventType = ((AliVHeader*)fInputEvent->GetHeader())->GetEventType();
381   if( fFiredTriggerClassName  !="" && !fAnaLED){
382     if(eventType!=7)
383       return kFALSE; //Only physics event, do not use for simulated events!!!
384     if(fDebug > 0) 
385       printf("AliCaloTrackReader::FillInputEvent() - FiredTriggerClass <%s>, selected class <%s>, compare name %d\n",
386              GetFiredTriggerClasses().Data(),fFiredTriggerClassName.Data(), GetFiredTriggerClasses().Contains(fFiredTriggerClassName));
387     if( !GetFiredTriggerClasses().Contains(fFiredTriggerClassName) ) return kFALSE;
388   }
389   else if(fAnaLED){
390 //        kStartOfRun =       1,    // START_OF_RUN
391 //        kEndOfRun =         2,    // END_OF_RUN
392 //        kStartOfRunFiles =  3,    // START_OF_RUN_FILES
393 //        kEndOfRunFiles =    4,    // END_OF_RUN_FILES
394 //        kStartOfBurst =     5,    // START_OF_BURST
395 //        kEndOfBurst =       6,    // END_OF_BURST
396 //        kPhysicsEvent =     7,    // PHYSICS_EVENT
397 //        kCalibrationEvent = 8,    // CALIBRATION_EVENT
398 //        kFormatError =      9,    // EVENT_FORMAT_ERROR
399 //        kStartOfData =      10,   // START_OF_DATA
400 //        kEndOfData =        11,   // END_OF_DATA
401 //        kSystemSoftwareTriggerEvent   = 12, // SYSTEM_SOFTWARE_TRIGGER_EVENT
402 //        kDetectorSoftwareTriggerEvent = 13  // DETECTOR_SOFTWARE_TRIGGER_EVENT
403          
404           if(eventType!=7 && fDebug > 1 )printf("AliCaloTrackReader::FillInputEvent() - DO LED, Event Type <%d>, 8 Calibration \n",  eventType);
405           if(eventType!=8)return kFALSE;
406   }
407                 
408   //In case of analysis of events with jets, skip those with jet pt > 5 pt hard 
409   if(fComparePtHardAndJetPt && GetStack()) {
410     if(!ComparePtHardAndJetPt()) return kFALSE ;
411   }
412
413   //In case of mixing events with other AOD file        
414  // if(fDataType == kAOD && fSecondInputAODTree){
415 //       
416 //    if(fDebug > 1) 
417 //      printf("AliCaloTrackReader::FillInputEvent() - Get event %d from second input AOD file \n", iEntry+fSecondInputFirstEvent);
418 //    if(fSecondInputAODTree->GetEntriesFast() <= iEntry+fSecondInputFirstEvent) {
419 //      if(fSecondInputAODTree->GetEntriesFast() == iEntry+fSecondInputFirstEvent) 
420 //                       printf("AliCaloTrackReader::FillInputEvent() - Skip events from event %d, no more events in second AOD file \n", iEntry);
421 //      return kFALSE;
422 //    }
423 //    
424 //    //Get the Event
425 //    Int_t nbytes = fSecondInputAODTree->GetEvent(iEntry+fSecondInputFirstEvent);
426 //    if ( nbytes == 0 ) {//If nothing in AOD
427 //      printf("AliCaloTrackReader::FillInputEvent() - Nothing in Second AOD input, STOP\n");
428 //      abort() ; 
429 //    }
430 //    
431 //  }
432         
433   //Fill Vertex array
434   
435   FillVertexArray();
436   
437   //Fill the arrays with cluster/tracks/cells data
438    if(fFillEMCALCells) 
439     FillInputEMCALCells();
440   if(fFillPHOSCells)  
441     FillInputPHOSCells();
442         
443   if(fFillCTS)   
444     FillInputCTS();
445   if(fFillEMCAL) 
446     FillInputEMCAL();
447   if(fFillPHOS)  
448     FillInputPHOS();
449
450   FillInputVZERO();
451         
452   return kTRUE ;
453 }
454
455 //__________________________________________________
456 void AliCaloTrackReader::ResetLists() {
457   //  Reset lists, called by the analysis maker 
458
459   if(fAODCTS)     fAODCTS     -> Clear();
460   if(fAODEMCAL)   fAODEMCAL   -> Clear();
461   if(fAODPHOS)    fAODPHOS    -> Clear();
462   if(fEMCALCells) fEMCALCells -> Clear();
463   if(fPHOSCells)  fPHOSCells  -> Clear();
464
465   fV0ADC[0] = 0;   fV0ADC[1] = 0; 
466   fV0Mul[0] = 0;   fV0Mul[1] = 0; 
467
468 }
469
470 //____________________________________________________________________________
471 void AliCaloTrackReader::SetInputEvent(AliVEvent* const input)  
472 {
473   fInputEvent  = input;
474   fMixedEvent = dynamic_cast<AliMixedEvent*>(GetInputEvent()) ; 
475   if (fMixedEvent) {
476     fNMixedEvent = fMixedEvent->GetNumberOfEvents() ; 
477   }
478
479   //Delete previous vertex
480   if(fVertex){
481     for (Int_t i = 0; i < fNMixedEvent; i++) {
482       delete [] fVertex[i] ; 
483     }
484     delete [] fVertex ;
485   }
486   
487   fVertex = new Double_t*[fNMixedEvent] ; 
488   for (Int_t i = 0; i < fNMixedEvent; i++) {
489     fVertex[i] = new Double_t[3] ; 
490     fVertex[i][0] = 0.0 ; 
491     fVertex[i][1] = 0.0 ; 
492     fVertex[i][2] = 0.0 ; 
493   }
494 }
495
496 //____________________________________________________________________________
497 void AliCaloTrackReader::GetVertex(Double_t vertex[3]) const {
498   //Return vertex position to be used for single event analysis
499   vertex[0]=fVertex[0][0];  
500   vertex[1]=fVertex[0][1];  
501   vertex[2]=fVertex[0][2];
502 }
503
504 //____________________________________________________________________________
505 void AliCaloTrackReader::GetVertex(Double_t vertex[3], const Int_t evtIndex) const {
506   //Return vertex position for mixed event, recover the vertex in a particular event.
507   
508   //Int_t evtIndex = 0; // for single events only one vertex stored in position 0, default value
509   //if (fMixedEvent && clusterID >=0) {
510   //  evtIndex=GetMixedEvent()->EventIndexForCaloCluster(clusterID) ; 
511   //}
512   
513   vertex[0]=fVertex[evtIndex][0];  vertex[1]=fVertex[evtIndex][1];  vertex[2]=fVertex[evtIndex][2];
514   
515 }
516 //
517
518
519 //____________________________________________________________________________
520 void AliCaloTrackReader::FillVertexArray() {
521   
522   //Fill data member with vertex
523   //In case of Mixed event, multiple vertices
524   
525   //Delete previous vertex
526   if(fVertex){
527     for (Int_t i = 0; i < fNMixedEvent; i++) {
528       delete [] fVertex[i] ; 
529     }
530     delete [] fVertex ;  
531   }
532   
533   fVertex = new Double_t*[fNMixedEvent] ; 
534   for (Int_t i = 0; i < fNMixedEvent; i++) {
535     fVertex[i] = new Double_t[3] ; 
536     fVertex[i][0] = 0.0 ; 
537     fVertex[i][1] = 0.0 ; 
538     fVertex[i][2] = 0.0 ; 
539   }          
540   
541   if (!fMixedEvent) { //Single event analysis
542     
543     if(fDataType!=kMC)fInputEvent->GetPrimaryVertex()->GetXYZ(fVertex[0]); 
544     else { 
545       fVertex[0][0]=0.;   fVertex[0][1]=0.;   fVertex[0][2]=0.;
546     }
547       
548     if(fDebug > 1)
549       printf("AliCaloTrackReader::FillVertexArray() - Single Event Vertex : %f,%f,%f\n",fVertex[0][0],fVertex[0][1],fVertex[0][2]);
550
551   } else { // MultiEvent analysis
552     for (Int_t iev = 0; iev < fNMixedEvent; iev++) {
553       if (fMixedEvent->GetVertexOfEvent(iev))
554         fMixedEvent->GetVertexOfEvent(iev)->GetXYZ(fVertex[iev]);
555       else { // no vertex found !!!!
556         AliWarning("No vertex found");
557       }
558
559       if(fDebug > 1)
560         printf("AliCaloTrackReader::FillVertexArray() - Multi Event %d Vertex : %f,%f,%f\n",iev,fVertex[iev][0],fVertex[iev][1],fVertex[iev][2]);
561
562     }
563   }
564   
565 }
566
567
568 //____________________________________________________________________________
569 void AliCaloTrackReader::FillInputCTS() {
570   //Return array with Central Tracking System (CTS) tracks
571   
572   if(fDebug > 2 ) printf("AliCaloTrackReader::FillInputCTS()\n");
573   
574   Int_t nTracks   = fInputEvent->GetNumberOfTracks() ;
575   Double_t p[3];
576   fTrackMult = 0;
577   Int_t nstatus = 0;
578   for (Int_t itrack =  0; itrack <  nTracks; itrack++) {////////////// track loop
579     AliVTrack * track = (AliVTrack*)fInputEvent->GetTrack(itrack) ; // retrieve track from esd
580
581     //Select tracks under certain conditions, TPCrefit, ITSrefit ... check the set bits
582     if (fTrackStatus && !((track->GetStatus() & fTrackStatus) == fTrackStatus)) 
583       continue ;
584     
585     nstatus++;
586     
587     if(fDataType==kESD && !fESDtrackCuts->AcceptTrack((AliESDtrack*)track)) continue;
588     
589     //Count the tracks in eta < 0.9
590     //printf("Eta %f cut  %f\n",TMath::Abs(track->Eta()),fTrackMultEtaCut);
591     if(TMath::Abs(track->Eta())< fTrackMultEtaCut) fTrackMult++;
592     
593     track->GetPxPyPz(p) ;
594     TLorentzVector momentum(p[0],p[1],p[2],0);
595     
596     if(fCTSPtMin < momentum.Pt()){
597       
598       if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"CTS")) 
599         continue;
600       
601       if(fDebug > 2 && momentum.Pt() > 0.1) 
602         printf("AliCaloTrackReader::FillInputCTS() - Selected tracks E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
603                momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
604       
605       if (fMixedEvent) {
606         track->SetID(itrack);
607       }
608       
609       fAODCTS->Add(track);        
610        
611     }//Pt and Fiducial cut passed. 
612   }// track loop
613         
614   //fAODCTSNormalInputEntries = fAODCTS->GetEntriesFast();
615   if(fDebug > 1) 
616     printf("AliCaloTrackReader::FillInputCTS()   - aod entries %d, input tracks %d, pass status %d, multipliticy %d\n", fAODCTS->GetEntriesFast(), nTracks, nstatus, fTrackMult);//fAODCTSNormalInputEntries);
617   
618     //  //If second input event available, add the clusters.
619     //  if(fSecondInputAODTree && fSecondInputAODEvent){
620     //    nTracks   = fSecondInputAODEvent->GetNumberOfTracks() ;
621     //    if(fDebug > 1) printf("AliCaloTrackReader::FillInputCTS()   - Add second input tracks, entries %d\n", nTracks) ;
622     //    for (Int_t itrack =  0; itrack <  nTracks; itrack++) {////////////// track loop
623     //            AliAODTrack * track = ((AliAODEvent*)fSecondInputAODEvent)->GetTrack(itrack) ; // retrieve track from esd
624     //            
625     //            //Select tracks under certain conditions, TPCrefit, ITSrefit ... check the set bits
626     //            if (fTrackStatus && !((track->GetStatus() & fTrackStatus) == fTrackStatus)) continue ;
627     //            
628     //            track->GetPxPyPz(p) ;
629     //            TLorentzVector momentum(p[0],p[1],p[2],0);
630     //            
631     //            if(fCTSPtMin < momentum.Pt()){
632     //
633     //                    if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"CTS")) continue;
634     //
635     //                    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",
636     //                                                                 momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
637     //                    
638     //                    fAODCTS->Add(track);
639     //                    
640     //            }//Pt and Fiducial cut passed. 
641     //    }// track loop
642     //    
643     //    if(fDebug > 1) printf("AliCaloTrackReader::FillInputCTS()   - aod normal entries %d, after second input %d\n", fAODCTSNormalInputEntries, fAODCTS->GetEntriesFast());
644     //  }       //second input loop
645     //  
646 }
647
648   //____________________________________________________________________________
649 void AliCaloTrackReader::FillInputEMCAL() {
650   //Return array with EMCAL clusters in aod format
651   
652   if(fDebug > 2 ) printf("AliCaloTrackReader::FillInputEMCAL()\n");
653   
654   //Loop to select clusters in fiducial cut and fill container with aodClusters
655   Int_t nclusters = fInputEvent->GetNumberOfCaloClusters();
656   for (Int_t iclus =  0; iclus <  nclusters; iclus++) {
657     AliVCluster * clus = 0;
658     if ( (clus = fInputEvent->GetCaloCluster(iclus)) ) {
659       if (IsEMCALCluster(clus)){
660         
661         //Check if the cluster contains any bad channel and if close to calorimeter borders
662         Int_t vindex = 0 ;  
663         if (fMixedEvent) 
664           vindex = fMixedEvent->EventIndexForCaloCluster(iclus);
665         
666         if(GetCaloUtils()->ClusterContainsBadChannel("EMCAL",clus->GetCellsAbsId(), clus->GetNCells())) 
667           continue;
668         if(!GetCaloUtils()->CheckCellFiducialRegion(clus, (AliVCaloCells*)fInputEvent->GetEMCALCells(), fInputEvent, vindex)) 
669           continue;
670         
671         TLorentzVector momentum ;
672         
673         clus->GetMomentum(momentum, fVertex[vindex]);      
674         
675         if(fEMCALPtMin < momentum.Pt()){
676           
677           if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"EMCAL")) 
678             continue;
679           
680           if(fDebug > 2 && momentum.E() > 0.1) 
681             printf("AliCaloTrackReader::FillInputEMCAL() - Selected clusters E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
682                    momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
683           
684           //Float_t pos[3];
685           //clus->GetPosition(pos);
686           //printf("Before Corrections: e %f, x %f, y %f, z %f\n",clus->E(),pos[0],pos[1],pos[2]);
687           
688           //Recalibrate the cluster energy 
689           if(GetCaloUtils()->IsRecalibrationOn()) {
690             Float_t energy = GetCaloUtils()->RecalibrateClusterEnergy(clus, GetEMCALCells());
691             clus->SetE(energy);
692             //printf("Recalibrated Energy %f\n",clus->E());  
693             GetCaloUtils()->RecalculateClusterShowerShapeParameters(GetEMCALCells(),clus);
694             GetCaloUtils()->RecalculateClusterPID(clus);
695
696           }
697           
698           //Correct non linearity
699           if(GetCaloUtils()->IsCorrectionOfClusterEnergyOn()){
700             GetCaloUtils()->CorrectClusterEnergy(clus) ;
701             //printf("Linearity Corrected Energy %f\n",clus->E());  
702           }
703           
704           //Recalculate cluster position
705           if(GetCaloUtils()->IsRecalculationOfClusterPositionOn()){
706             GetCaloUtils()->RecalculateClusterPosition(GetEMCALCells(),clus); 
707             //clus->GetPosition(pos);
708             //printf("After  Corrections: e %f, x %f, y %f, z %f\n",clus->E(),pos[0],pos[1],pos[2]);
709           }
710           
711           if (fMixedEvent) {
712             clus->SetID(iclus) ; 
713           }
714             
715           fAODEMCAL->Add(clus); 
716           
717         }//Pt and Fiducial cut passed.
718       }//EMCAL cluster
719     }// cluster exists
720   }// cluster loop
721   
722   //Recalculate track matching
723   GetCaloUtils()->RecalculateClusterTrackMatching(fInputEvent);
724   
725   //fAODEMCALNormalInputEntries = fAODEMCAL->GetEntriesFast();
726   if(fDebug > 1) printf("AliCaloTrackReader::FillInputEMCAL() - aod entries %d\n",  fAODEMCAL->GetEntriesFast());//fAODEMCALNormalInputEntries);
727   
728     //If second input event available, add the clusters.
729     //  if(fSecondInputAODTree && fSecondInputAODEvent){
730     //    GetSecondInputAODVertex(v);
731     //    nclusters = ((AliAODEvent*)fSecondInputAODEvent)->GetNumberOfCaloClusters();
732     //    if(fDebug > 1) printf("AliCaloTrackReader::FillInputEMCAL() - Add second input clusters, entries %d\n", nclusters) ;
733     //          for (Int_t iclus =  0; iclus < nclusters; iclus++) {
734     //                  AliAODCaloCluster * clus = 0;
735     //                  if ( (clus = ((AliAODEvent*)fSecondInputAODEvent)->GetCaloCluster(iclus)) ) {
736     //                          if (clus->IsEMCAL()){
737     //                                  TLorentzVector momentum ;
738     //                                  clus->GetMomentum(momentum, v);      
739     //                                  
740     //                                  if(fEMCALPtMin < momentum.Pt()){
741     //
742     //                                          if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"EMCAL")) continue;
743     //
744     //                                          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",
745     //                                                                                                                                  momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
746     //                                    fAODEMCAL->Add(clus); 
747     //                                  }//Pt and Fiducial cut passed.
748     //                          }//EMCAL cluster
749     //                  }// cluster exists
750     //          }// cluster loop
751     //          
752     //    if(fDebug > 1) printf("AliCaloTrackReader::FillInputEMCAL() - aod normal entries %d, after second input %d\n", fAODEMCALNormalInputEntries, fAODEMCAL->GetEntriesFast());
753     //
754     //  } //second input loop
755 }
756
757   //____________________________________________________________________________
758 void AliCaloTrackReader::FillInputPHOS() {
759   //Return array with PHOS clusters in aod format
760   
761   if(fDebug > 2 ) printf("AliCaloTrackReader::FillInputPHOS()\n");
762           
763   //Loop to select clusters in fiducial cut and fill container with aodClusters
764   Int_t nclusters = fInputEvent->GetNumberOfCaloClusters();
765   for (Int_t iclus = 0; iclus < nclusters; iclus++) {
766     AliVCluster * clus = 0;
767     if ( (clus = fInputEvent->GetCaloCluster(iclus)) ) {
768       if (IsPHOSCluster(clus)){
769         //Check if the cluster contains any bad channel and if close to calorimeter borders
770         Int_t vindex = 0 ;  
771         if (fMixedEvent) 
772           vindex = fMixedEvent->EventIndexForCaloCluster(iclus);
773         if( GetCaloUtils()->ClusterContainsBadChannel("PHOS",clus->GetCellsAbsId(), clus->GetNCells())) 
774           continue;
775         if(!GetCaloUtils()->CheckCellFiducialRegion(clus, fInputEvent->GetPHOSCells(), fInputEvent, vindex)) 
776           continue;
777         
778         TLorentzVector momentum ;
779         
780         clus->GetMomentum(momentum, fVertex[vindex]);      
781         
782         if(fPHOSPtMin < momentum.Pt()){
783           
784           if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"PHOS")) 
785             continue;
786           
787           if(fDebug > 2 && momentum.E() > 0.1) 
788             printf("AliCaloTrackReader::FillInputPHOS() - Selected clusters E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
789                    momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
790           
791             //Recalibrate the cluster energy 
792           if(GetCaloUtils()->IsRecalibrationOn()) {
793             Float_t energy = GetCaloUtils()->RecalibrateClusterEnergy(clus, (AliAODCaloCells*)GetPHOSCells());
794             clus->SetE(energy);
795           }
796           
797           if (fMixedEvent) {
798             clus->SetID(iclus) ; 
799           }              
800           
801           fAODPHOS->Add(clus);  
802           
803         }//Pt and Fiducial cut passed.
804       }//PHOS cluster
805     }//cluster exists
806   }//esd cluster loop
807   
808   //fAODPHOSNormalInputEntries = fAODPHOS->GetEntriesFast() ;
809   if(fDebug > 1) printf("AliCaloTrackReader::FillInputPHOS()  - aod entries %d\n",  fAODPHOS->GetEntriesFast());//fAODPHOSNormalInputEntries);
810   
811     //If second input event available, add the clusters.
812     //  if(fSecondInputAODTree && fSecondInputAODEvent){  
813     //    GetSecondInputAODVertex(v);
814     //    nclusters = ((AliAODEvent*)fSecondInputAODEvent)->GetNumberOfCaloClusters();
815     //    if(fDebug > 1) printf("AliCaloTrackReader::FillInputPHOS()  - Add second input clusters, entries %d\n", nclusters);
816     //          for (Int_t iclus =  0; iclus < nclusters; iclus++) {
817     //                  AliAODCaloCluster * clus = 0;
818     //                  if ( (clus = ((AliAODEvent*)fSecondInputAODEvent)->GetCaloCluster(iclus)) ) {
819     //                          if (clus->IsPHOS()){
820     //                                  TLorentzVector momentum ;
821     //                                  clus->GetMomentum(momentum, v);      
822     //                                  
823     //                                  if(fPHOSPtMin < momentum.Pt()){
824     //
825     //                                          if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"PHOS")) continue;
826     //
827     //                                          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",
828     //                                                                                                                                  momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
829     //                                          fAODPHOS->Add(clus);    
830     //                                  }//Pt and Fiducial cut passed.
831     //                          }//PHOS cluster
832     //                  }// cluster exists
833     //          }// cluster loop
834     //          if(fDebug > 1) printf("AliCaloTrackReader::FillInputPHOS()  - aod normal entries %d, after second input %d\n", fAODPHOSNormalInputEntries, fAODPHOS->GetEntriesFast());
835     //  }       //second input loop
836   
837 }
838
839 //____________________________________________________________________________
840 void AliCaloTrackReader::FillInputEMCALCells() {
841     //Return array with EMCAL cells in aod format
842   
843   fEMCALCells = fInputEvent->GetEMCALCells(); 
844   
845 }
846
847 //____________________________________________________________________________
848 void AliCaloTrackReader::FillInputPHOSCells() {
849     //Return array with PHOS cells in aod format
850   
851   fPHOSCells = fInputEvent->GetPHOSCells(); 
852   
853 }
854
855
856 //____________________________________________________________________________
857 Bool_t AliCaloTrackReader::IsEMCALCluster(AliVCluster* cluster) const {
858   // Check if it is a cluster from EMCAL. For old AODs cluster type has
859   // different number and need to patch here
860     
861   if(fDataType==kAOD && fOldAOD)
862   {
863     if (cluster->GetType() == 2) return kTRUE;
864     else                         return kFALSE;
865   }
866   else 
867   {
868     return cluster->IsEMCAL();
869   }
870
871 }
872
873 //____________________________________________________________________________
874 Bool_t AliCaloTrackReader::IsPHOSCluster(AliVCluster * cluster) const {
875   //Check if it is a cluster from PHOS.For old AODs cluster type has
876   // different number and need to patch here
877   
878   if(fDataType==kAOD && fOldAOD)
879   {
880     Int_t type = cluster->GetType();
881     if (type == 0 || type == 1) return kTRUE;
882     else                        return kFALSE;
883   }
884   else 
885   {
886     return cluster->IsPHOS();
887   }
888   
889 }
890