]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALRecoUtils.cxx
Merge remote-tracking branch 'origin/master' into mergingFlat
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALRecoUtils.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
16 /* $Id: AliEMCALRecoUtils.cxx | Sun Dec 8 06:56:48 2013 +0100 | Constantin Loizides  $ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //
20 // Class AliEMCALRecoUtils
21 // Some utilities to recalculate the cluster position or energy linearity
22 //
23 //
24 // Author:  Gustavo Conesa (LPSC- Grenoble) 
25 //          Track matching part: Rongrong Ma (Yale)
26
27 ///////////////////////////////////////////////////////////////////////////////
28 // --- standard c ---
29
30 // standard C++ includes
31 //#include <Riostream.h>
32
33 // ROOT includes
34 #include <TGeoManager.h>
35 #include <TGeoMatrix.h>
36 #include <TGeoBBox.h>
37 #include <TH2F.h>
38 #include <TArrayI.h>
39 #include <TArrayF.h>
40 #include <TObjArray.h>
41
42 // STEER includes
43 #include "AliVCluster.h"
44 #include "AliVCaloCells.h"
45 #include "AliLog.h"
46 #include "AliPID.h"
47 #include "AliESDEvent.h"
48 #include "AliAODEvent.h"
49 #include "AliESDtrack.h"
50 #include "AliAODTrack.h"
51 #include "AliExternalTrackParam.h"
52 #include "AliESDfriendTrack.h"
53 #include "AliTrackerBase.h"
54
55 // EMCAL includes
56 #include "AliEMCALRecoUtils.h"
57 #include "AliEMCALGeometry.h"
58 #include "AliTrackerBase.h"
59 #include "AliEMCALPIDUtils.h"
60
61 ClassImp(AliEMCALRecoUtils)
62   
63 //_____________________________________
64 AliEMCALRecoUtils::AliEMCALRecoUtils():
65   fParticleType(0),                       fPosAlgo(0),                            fW0(0), 
66   fNonLinearityFunction(0),               fNonLinearThreshold(0),
67   fSmearClusterEnergy(kFALSE),            fRandom(),
68   fCellsRecalibrated(kFALSE),             fRecalibration(kFALSE),                 fEMCALRecalibrationFactors(),
69   fTimeRecalibration(kFALSE),             fEMCALTimeRecalibrationFactors(),       fUseRunCorrectionFactors(kFALSE),       
70   fRemoveBadChannels(kFALSE),             fRecalDistToBadChannels(kFALSE),        fEMCALBadChannelMap(),
71   fNCellsFromEMCALBorder(0),              fNoEMCALBorderAtEta0(kTRUE),
72   fRejectExoticCluster(kFALSE),           fRejectExoticCells(kFALSE), 
73   fExoticCellFraction(0),                 fExoticCellDiffTime(0),                 fExoticCellMinAmplitude(0),
74   fPIDUtils(),                            fAODFilterMask(0),
75   fAODHybridTracks(0),                    fAODTPCOnlyTracks(0),
76   fMatchedTrackIndex(0x0),                fMatchedClusterIndex(0x0), 
77   fResidualEta(0x0), fResidualPhi(0x0),   fCutEtaPhiSum(kFALSE),                  fCutEtaPhiSeparate(kFALSE), 
78   fCutR(0),                               fCutEta(0),                             fCutPhi(0),
79   fClusterWindow(0),                      fMass(0),                           
80   fStepSurface(0),                        fStepCluster(0),
81   fITSTrackSA(kFALSE),                    fEMCalSurfaceDistance(440.),
82   fTrackCutsType(0),                      fCutMinTrackPt(0),                      fCutMinNClusterTPC(0), 
83   fCutMinNClusterITS(0),                  fCutMaxChi2PerClusterTPC(0),            fCutMaxChi2PerClusterITS(0),
84   fCutRequireTPCRefit(kFALSE),            fCutRequireITSRefit(kFALSE),            fCutAcceptKinkDaughters(kFALSE),
85   fCutMaxDCAToVertexXY(0),                fCutMaxDCAToVertexZ(0),                 fCutDCAToVertex2D(kFALSE),
86   fCutRequireITSStandAlone(kFALSE),       fCutRequireITSpureSA(kFALSE)
87 {
88 //
89   // Constructor.
90   // Initialize all constant values which have to be used
91   // during Reco algorithm execution
92   //
93   
94   // Init parameters
95   InitParameters();
96   
97   //Track matching
98   fMatchedTrackIndex     = new TArrayI();
99   fMatchedClusterIndex   = new TArrayI();
100   fResidualPhi           = new TArrayF();
101   fResidualEta           = new TArrayF();
102   fPIDUtils              = new AliEMCALPIDUtils();
103
104 }
105
106 //______________________________________________________________________
107 AliEMCALRecoUtils::AliEMCALRecoUtils(const AliEMCALRecoUtils & reco) 
108 : TNamed(reco), 
109   fParticleType(reco.fParticleType),                         fPosAlgo(reco.fPosAlgo),     fW0(reco.fW0),
110   fNonLinearityFunction(reco.fNonLinearityFunction),         fNonLinearThreshold(reco.fNonLinearThreshold),
111   fSmearClusterEnergy(reco.fSmearClusterEnergy),             fRandom(),
112   fCellsRecalibrated(reco.fCellsRecalibrated),
113   fRecalibration(reco.fRecalibration),                       fEMCALRecalibrationFactors(reco.fEMCALRecalibrationFactors),
114   fTimeRecalibration(reco.fTimeRecalibration),               fEMCALTimeRecalibrationFactors(reco.fEMCALTimeRecalibrationFactors),
115   fUseRunCorrectionFactors(reco.fUseRunCorrectionFactors),   
116   fRemoveBadChannels(reco.fRemoveBadChannels),               fRecalDistToBadChannels(reco.fRecalDistToBadChannels),
117   fEMCALBadChannelMap(reco.fEMCALBadChannelMap),
118   fNCellsFromEMCALBorder(reco.fNCellsFromEMCALBorder),       fNoEMCALBorderAtEta0(reco.fNoEMCALBorderAtEta0),
119   fRejectExoticCluster(reco.fRejectExoticCluster),           fRejectExoticCells(reco.fRejectExoticCells), 
120   fExoticCellFraction(reco.fExoticCellFraction),             fExoticCellDiffTime(reco.fExoticCellDiffTime),               
121   fExoticCellMinAmplitude(reco.fExoticCellMinAmplitude),
122   fPIDUtils(reco.fPIDUtils),                                 fAODFilterMask(reco.fAODFilterMask),
123   fAODHybridTracks(reco.fAODHybridTracks),                   fAODTPCOnlyTracks(reco.fAODTPCOnlyTracks),
124   fMatchedTrackIndex(  reco.fMatchedTrackIndex?  new TArrayI(*reco.fMatchedTrackIndex):0x0),
125   fMatchedClusterIndex(reco.fMatchedClusterIndex?new TArrayI(*reco.fMatchedClusterIndex):0x0),
126   fResidualEta(        reco.fResidualEta?        new TArrayF(*reco.fResidualEta):0x0),
127   fResidualPhi(        reco.fResidualPhi?        new TArrayF(*reco.fResidualPhi):0x0),
128   fCutEtaPhiSum(reco.fCutEtaPhiSum),                         fCutEtaPhiSeparate(reco.fCutEtaPhiSeparate), 
129   fCutR(reco.fCutR),        fCutEta(reco.fCutEta),           fCutPhi(reco.fCutPhi),
130   fClusterWindow(reco.fClusterWindow),
131   fMass(reco.fMass),        fStepSurface(reco.fStepSurface), fStepCluster(reco.fStepCluster),
132   fITSTrackSA(reco.fITSTrackSA),                             fEMCalSurfaceDistance(440.),
133   fTrackCutsType(reco.fTrackCutsType),                       fCutMinTrackPt(reco.fCutMinTrackPt), 
134   fCutMinNClusterTPC(reco.fCutMinNClusterTPC),               fCutMinNClusterITS(reco.fCutMinNClusterITS), 
135   fCutMaxChi2PerClusterTPC(reco.fCutMaxChi2PerClusterTPC),   fCutMaxChi2PerClusterITS(reco.fCutMaxChi2PerClusterITS),
136   fCutRequireTPCRefit(reco.fCutRequireTPCRefit),             fCutRequireITSRefit(reco.fCutRequireITSRefit),
137   fCutAcceptKinkDaughters(reco.fCutAcceptKinkDaughters),     fCutMaxDCAToVertexXY(reco.fCutMaxDCAToVertexXY),    
138   fCutMaxDCAToVertexZ(reco.fCutMaxDCAToVertexZ),             fCutDCAToVertex2D(reco.fCutDCAToVertex2D),
139   fCutRequireITSStandAlone(reco.fCutRequireITSStandAlone),   fCutRequireITSpureSA(reco.fCutRequireITSpureSA)
140 {
141   //Copy ctor
142   
143   for (Int_t i = 0; i < 15 ; i++) { fMisalRotShift[i]      = reco.fMisalRotShift[i]      ; 
144                                    fMisalTransShift[i]    = reco.fMisalTransShift[i]    ; }
145   for (Int_t i = 0; i < 7  ; i++) { fNonLinearityParams[i] = reco.fNonLinearityParams[i] ; }
146   for (Int_t i = 0; i < 3  ; i++) { fSmearClusterParam[i]  = reco.fSmearClusterParam[i]  ; }
147
148 }
149
150
151 //______________________________________________________________________
152 AliEMCALRecoUtils & AliEMCALRecoUtils::operator = (const AliEMCALRecoUtils & reco) 
153 {
154   //Assignment operator
155   
156   if (this == &reco)return *this;
157   ((TNamed *)this)->operator=(reco);
158
159   for (Int_t i = 0; i < 15 ; i++) { fMisalTransShift[i]    = reco.fMisalTransShift[i]    ; 
160                                    fMisalRotShift[i]      = reco.fMisalRotShift[i]      ; }
161   for (Int_t i = 0; i < 7  ; i++) { fNonLinearityParams[i] = reco.fNonLinearityParams[i] ; }
162   for (Int_t i = 0; i < 3  ; i++) { fSmearClusterParam[i]  = reco.fSmearClusterParam[i]  ; }   
163   
164   fParticleType              = reco.fParticleType;
165   fPosAlgo                   = reco.fPosAlgo; 
166   fW0                        = reco.fW0;
167   
168   fNonLinearityFunction      = reco.fNonLinearityFunction;
169   fNonLinearThreshold        = reco.fNonLinearThreshold;
170   fSmearClusterEnergy        = reco.fSmearClusterEnergy;
171
172   fCellsRecalibrated         = reco.fCellsRecalibrated;
173   fRecalibration             = reco.fRecalibration;
174   fEMCALRecalibrationFactors = reco.fEMCALRecalibrationFactors;
175
176   fTimeRecalibration             = reco.fTimeRecalibration;
177   fEMCALTimeRecalibrationFactors = reco.fEMCALTimeRecalibrationFactors;
178
179   fUseRunCorrectionFactors   = reco.fUseRunCorrectionFactors;
180   
181   fRemoveBadChannels         = reco.fRemoveBadChannels;
182   fRecalDistToBadChannels    = reco.fRecalDistToBadChannels;
183   fEMCALBadChannelMap        = reco.fEMCALBadChannelMap;
184   
185   fNCellsFromEMCALBorder     = reco.fNCellsFromEMCALBorder;
186   fNoEMCALBorderAtEta0       = reco.fNoEMCALBorderAtEta0;
187   
188   fRejectExoticCluster       = reco.fRejectExoticCluster;           
189   fRejectExoticCells         = reco.fRejectExoticCells; 
190   fExoticCellFraction        = reco.fExoticCellFraction;
191   fExoticCellDiffTime        = reco.fExoticCellDiffTime;              
192   fExoticCellMinAmplitude    = reco.fExoticCellMinAmplitude;
193   
194   fPIDUtils                  = reco.fPIDUtils;
195
196   fAODFilterMask             = reco.fAODFilterMask;
197   fAODHybridTracks           = reco.fAODHybridTracks;
198   fAODTPCOnlyTracks          = reco.fAODTPCOnlyTracks;
199   
200   fCutEtaPhiSum              = reco.fCutEtaPhiSum;
201   fCutEtaPhiSeparate         = reco.fCutEtaPhiSeparate;
202   fCutR                      = reco.fCutR;
203   fCutEta                    = reco.fCutEta;
204   fCutPhi                    = reco.fCutPhi;
205   fClusterWindow             = reco.fClusterWindow;
206   fMass                      = reco.fMass;
207   fStepSurface               = reco.fStepSurface;
208   fStepCluster               = reco.fStepCluster;
209   fITSTrackSA                = reco.fITSTrackSA;
210   fEMCalSurfaceDistance      = reco.fEMCalSurfaceDistance;
211   
212   fTrackCutsType             = reco.fTrackCutsType;
213   fCutMinTrackPt             = reco.fCutMinTrackPt;
214   fCutMinNClusterTPC         = reco.fCutMinNClusterTPC;
215   fCutMinNClusterITS         = reco.fCutMinNClusterITS; 
216   fCutMaxChi2PerClusterTPC   = reco.fCutMaxChi2PerClusterTPC;
217   fCutMaxChi2PerClusterITS   = reco.fCutMaxChi2PerClusterITS;
218   fCutRequireTPCRefit        = reco.fCutRequireTPCRefit;
219   fCutRequireITSRefit        = reco.fCutRequireITSRefit;
220   fCutAcceptKinkDaughters    = reco.fCutAcceptKinkDaughters;
221   fCutMaxDCAToVertexXY       = reco.fCutMaxDCAToVertexXY;
222   fCutMaxDCAToVertexZ        = reco.fCutMaxDCAToVertexZ;
223   fCutDCAToVertex2D          = reco.fCutDCAToVertex2D;
224   fCutRequireITSStandAlone   = reco.fCutRequireITSStandAlone; 
225   fCutRequireITSpureSA       = reco.fCutRequireITSpureSA;
226
227   if (reco.fResidualEta) {
228     // assign or copy construct
229     if (fResidualEta) { 
230       *fResidualEta = *reco.fResidualEta;
231     } else {
232       fResidualEta = new TArrayF(*reco.fResidualEta);
233     }
234   } else {
235     if (fResidualEta) delete fResidualEta;
236     fResidualEta = 0;
237   }
238   
239   if (reco.fResidualPhi) {
240     // assign or copy construct
241     if (fResidualPhi) { 
242       *fResidualPhi = *reco.fResidualPhi;
243     } else {
244       fResidualPhi = new TArrayF(*reco.fResidualPhi);
245     }
246   } else {
247     if (fResidualPhi) delete fResidualPhi;
248     fResidualPhi = 0;
249   }
250   
251   if (reco.fMatchedTrackIndex) {
252     // assign or copy construct
253     if (fMatchedTrackIndex) { 
254       *fMatchedTrackIndex = *reco.fMatchedTrackIndex;
255     } else  { 
256       fMatchedTrackIndex = new TArrayI(*reco.fMatchedTrackIndex);
257     }
258   } else {
259     if (fMatchedTrackIndex) delete fMatchedTrackIndex;
260     fMatchedTrackIndex = 0;
261   }  
262   
263   if (reco.fMatchedClusterIndex) {
264     // assign or copy construct
265     if (fMatchedClusterIndex) { 
266       *fMatchedClusterIndex = *reco.fMatchedClusterIndex;
267     } else {
268       fMatchedClusterIndex = new TArrayI(*reco.fMatchedClusterIndex);
269     }
270   } else {
271     if (fMatchedClusterIndex) delete fMatchedClusterIndex;
272     fMatchedClusterIndex = 0;
273   }
274   
275   return *this;
276 }
277
278 //_____________________________________
279 AliEMCALRecoUtils::~AliEMCALRecoUtils()
280 {
281   //Destructor.
282   
283   if (fEMCALRecalibrationFactors) { 
284     fEMCALRecalibrationFactors->Clear();
285     delete fEMCALRecalibrationFactors;
286   }  
287   
288   if (fEMCALTimeRecalibrationFactors) { 
289     fEMCALTimeRecalibrationFactors->Clear();
290     delete fEMCALTimeRecalibrationFactors;
291   }  
292   
293   if (fEMCALBadChannelMap) { 
294     fEMCALBadChannelMap->Clear();
295     delete fEMCALBadChannelMap;
296   }
297  
298   delete fMatchedTrackIndex   ; 
299   delete fMatchedClusterIndex ; 
300   delete fResidualEta         ; 
301   delete fResidualPhi         ; 
302   delete fPIDUtils            ;
303
304   InitTrackCuts();
305 }
306
307 //_______________________________________________________________________________
308 Bool_t AliEMCALRecoUtils::AcceptCalibrateCell(Int_t absID, Int_t bc,
309                                               Float_t  & amp,    Double_t & time, 
310                                               AliVCaloCells* cells) 
311 {
312   // Reject cell if criteria not passed and calibrate it
313   
314   AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
315   
316   if (absID < 0 || absID >= 24*48*geom->GetNumberOfSuperModules()) 
317     return kFALSE;
318   
319   Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1; 
320   
321   if (!geom->GetCellIndex(absID,imod,iTower,iIphi,iIeta)) {
322     // cell absID does not exist
323     amp=0; time = 1.e9;
324     return kFALSE; 
325   }
326   
327   geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);  
328
329   // Do not include bad channels found in analysis,
330   if (IsBadChannelsRemovalSwitchedOn() && GetEMCALChannelStatus(imod, ieta, iphi)) {
331     return kFALSE;
332   }
333   
334   //Recalibrate energy
335   amp  = cells->GetCellAmplitude(absID);
336   if (!fCellsRecalibrated && IsRecalibrationOn())
337     amp *= GetEMCALChannelRecalibrationFactor(imod,ieta,iphi);
338   
339   // Recalibrate time
340   time = cells->GetCellTime(absID);
341   
342   RecalibrateCellTime(absID,bc,time);
343   
344   return kTRUE;
345 }
346
347 //_____________________________________________________________________________
348 Bool_t AliEMCALRecoUtils::CheckCellFiducialRegion(const AliEMCALGeometry* geom, 
349                                                   const AliVCluster* cluster, 
350                                                   AliVCaloCells* cells) 
351 {
352   // Given the list of AbsId of the cluster, get the maximum cell and 
353   // check if there are fNCellsFromBorder from the calorimeter border
354   
355   if (!cluster)
356   {
357     AliInfo("Cluster pointer null!");
358     return kFALSE;
359   }
360   
361   //If the distance to the border is 0 or negative just exit accept all clusters
362   if (cells->GetType()==AliVCaloCells::kEMCALCell && fNCellsFromEMCALBorder <= 0 ) 
363     return kTRUE;
364   
365   Int_t absIdMax  = -1, iSM =-1, ieta = -1, iphi = -1;
366   Bool_t shared = kFALSE;
367   GetMaxEnergyCell(geom, cells, cluster, absIdMax,  iSM, ieta, iphi, shared);
368   
369   AliDebug(2,Form("Cluster Max AbsId %d, Cell Energy %2.2f, Cluster Energy %2.2f, Ncells from border %d, EMCAL eta=0 %d\n", 
370                   absIdMax, cells->GetCellAmplitude(absIdMax), cluster->E(), fNCellsFromEMCALBorder, fNoEMCALBorderAtEta0));
371   
372   if (absIdMax==-1) return kFALSE;
373   
374   //Check if the cell is close to the borders:
375   Bool_t okrow = kFALSE;
376   Bool_t okcol = kFALSE;
377   
378   if (iSM < 0 || iphi < 0 || ieta < 0 ) {
379     AliFatal(Form("Negative value for super module: %d, or cell ieta: %d, or cell iphi: %d, check EMCAL geometry name\n",
380                   iSM,ieta,iphi));
381     return kFALSE; // trick coverity
382   }
383   
384   //Check rows/phi
385   Int_t iPhiLast = 24;
386    if( geom->GetSMType(iSM) == AliEMCALGeometry::kEMCAL_Half ) iPhiLast /= 2;
387    else if (  geom->GetSMType(iSM) == AliEMCALGeometry::kEMCAL_3rd ) iPhiLast /= 3;// 1/3 sm case
388   
389   if(iphi >= fNCellsFromEMCALBorder && iphi < iPhiLast - fNCellsFromEMCALBorder) okrow = kTRUE; 
390
391   //Check columns/eta
392   Int_t iEtaLast = 48;
393   if(!fNoEMCALBorderAtEta0 || geom->IsDCALSM(iSM)) {// conside inner border
394      if(  geom->GetSMType(iSM) == AliEMCALGeometry::kDCAL_Standard )  iEtaLast = iEtaLast*2/3;        
395      if(ieta  > fNCellsFromEMCALBorder && ieta < iEtaLast-fNCellsFromEMCALBorder) okcol = kTRUE;  
396   } else {
397     if (iSM%2==0) {
398      if (ieta >= fNCellsFromEMCALBorder)     okcol = kTRUE;  
399     } else {
400      if(ieta <  iEtaLast-fNCellsFromEMCALBorder)  okcol = kTRUE; 
401     }
402   }//eta 0 not checked
403   
404   AliDebug(2,Form("EMCAL Cluster in %d cells fiducial volume: ieta %d, iphi %d, SM %d:  column? %d, row? %d\nq",
405                   fNCellsFromEMCALBorder, ieta, iphi, iSM, okcol, okrow));
406   
407   if (okcol && okrow) {
408     //printf("Accept\n");
409     return kTRUE;
410   } else  {
411     //printf("Reject\n");
412     AliDebug(2,Form("Reject cluster in border, max cell : ieta %d, iphi %d, SM %d\n",ieta, iphi, iSM));
413     return kFALSE;
414   }
415 }  
416
417 //_______________________________________________________________________________
418 Bool_t AliEMCALRecoUtils::ClusterContainsBadChannel(const AliEMCALGeometry* geom, 
419                                                     const UShort_t* cellList, 
420                                                     Int_t nCells)
421 {
422   // Check that in the cluster cells, there is no bad channel of those stored 
423   // in fEMCALBadChannelMap or fPHOSBadChannelMap
424   
425   if (!fRemoveBadChannels)  return kFALSE;
426   if (!fEMCALBadChannelMap) return kFALSE;
427   
428   Int_t icol = -1;
429   Int_t irow = -1;
430   Int_t imod = -1;
431   for (Int_t iCell = 0; iCell<nCells; iCell++) {
432     //Get the column and row
433     Int_t iTower = -1, iIphi = -1, iIeta = -1; 
434     geom->GetCellIndex(cellList[iCell],imod,iTower,iIphi,iIeta); 
435     if (fEMCALBadChannelMap->GetEntries() <= imod) continue;
436     geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,irow,icol);      
437     if (GetEMCALChannelStatus(imod, icol, irow)) {
438       AliDebug(2,Form("Cluster with bad channel: SM %d, col %d, row %d\n",imod, icol, irow));
439       return kTRUE;
440     }
441   }// cell cluster loop
442
443   return kFALSE;
444 }
445
446
447 //___________________________________________________________________________
448 Float_t AliEMCALRecoUtils::GetECross(Int_t absID, Double_t tcell,
449                                      AliVCaloCells* cells, Int_t bc)
450 {
451   //Calculate the energy in the cross around the energy given cell
452   
453   AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
454   
455   Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1; 
456   geom->GetCellIndex(absID,imod,iTower,iIphi,iIeta); 
457   geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);  
458   
459   //Get close cells index, energy and time, not in corners
460   
461   Int_t absID1 = -1;
462   Int_t absID2 = -1;
463   
464   if ( iphi < AliEMCALGeoParams::fgkEMCALRows-1) absID1 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi+1, ieta);
465   if ( iphi > 0 )                                absID2 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi-1, ieta);
466   
467   // In case of cell in eta = 0 border, depending on SM shift the cross cell index
468   
469   Int_t absID3 = -1;
470   Int_t absID4 = -1;
471   
472   if ( ieta == AliEMCALGeoParams::fgkEMCALCols-1 && !(imod%2) ) {
473     absID3 = geom-> GetAbsCellIdFromCellIndexes(imod+1, iphi, 0);
474     absID4 = geom-> GetAbsCellIdFromCellIndexes(imod,   iphi, ieta-1); 
475   } else if ( ieta == 0 && imod%2 ) {
476     absID3 = geom-> GetAbsCellIdFromCellIndexes(imod,   iphi, ieta+1);
477     absID4 = geom-> GetAbsCellIdFromCellIndexes(imod-1, iphi, AliEMCALGeoParams::fgkEMCALCols-1); 
478   } else {
479     if ( ieta < AliEMCALGeoParams::fgkEMCALCols-1 ) 
480       absID3 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta+1);
481     if ( ieta > 0 )                                 
482       absID4 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta-1); 
483   }
484   
485   //printf("IMOD %d, AbsId %d, a %d, b %d, c %d e %d \n",imod,absID,absID1,absID2,absID3,absID4);
486   
487   Float_t  ecell1  = 0, ecell2  = 0, ecell3  = 0, ecell4  = 0;
488   Double_t tcell1  = 0, tcell2  = 0, tcell3  = 0, tcell4  = 0;
489
490   AcceptCalibrateCell(absID1,bc, ecell1,tcell1,cells); 
491   AcceptCalibrateCell(absID2,bc, ecell2,tcell2,cells); 
492   AcceptCalibrateCell(absID3,bc, ecell3,tcell3,cells); 
493   AcceptCalibrateCell(absID4,bc, ecell4,tcell4,cells); 
494   
495   if (TMath::Abs(tcell-tcell1)*1.e9 > fExoticCellDiffTime) ecell1 = 0 ;
496   if (TMath::Abs(tcell-tcell2)*1.e9 > fExoticCellDiffTime) ecell2 = 0 ;
497   if (TMath::Abs(tcell-tcell3)*1.e9 > fExoticCellDiffTime) ecell3 = 0 ;
498   if (TMath::Abs(tcell-tcell4)*1.e9 > fExoticCellDiffTime) ecell4 = 0 ;
499   
500   return ecell1+ecell2+ecell3+ecell4;
501 }
502
503 //_____________________________________________________________________________________________
504 Bool_t AliEMCALRecoUtils::IsExoticCell(Int_t absID, AliVCaloCells* cells, Int_t bc)
505 {
506   // Look to cell neighbourhood and reject if it seems exotic
507   // Do before recalibrating the cells
508
509   if (!fRejectExoticCells) return kFALSE;
510   
511   Float_t  ecell  = 0;
512   Double_t tcell  = 0;
513   Bool_t   accept = AcceptCalibrateCell(absID, bc, ecell ,tcell ,cells); 
514   
515   if (!accept) return kTRUE; // reject this cell
516   
517   if (ecell < fExoticCellMinAmplitude) return kFALSE; // do not reject low energy cells
518
519   Float_t eCross = GetECross(absID,tcell,cells,bc);
520   
521   if (1-eCross/ecell > fExoticCellFraction) {
522     AliDebug(2,Form("AliEMCALRecoUtils::IsExoticCell() - EXOTIC CELL id %d, eCell %f, eCross %f, 1-eCross/eCell %f\n",
523                     absID,ecell,eCross,1-eCross/ecell));
524     return kTRUE;
525   }
526
527   return kFALSE;
528 }
529
530 //___________________________________________________________________
531 Bool_t AliEMCALRecoUtils::IsExoticCluster(const AliVCluster *cluster, 
532                                           AliVCaloCells *cells, 
533                                           Int_t bc) 
534 {
535   // Check if the cluster highest energy tower is exotic
536   
537   if (!cluster) {
538     AliInfo("Cluster pointer null!");
539     return kFALSE;
540   }
541   
542   if (!fRejectExoticCluster) return kFALSE;
543   
544   // Get highest energy tower
545   AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
546   Int_t iSupMod = -1, absId = -1, ieta = -1, iphi = -1;
547   Bool_t shared = kFALSE;
548   GetMaxEnergyCell(geom, cells, cluster, absId, iSupMod, ieta, iphi, shared);
549   
550   return IsExoticCell(absId,cells,bc);
551 }
552
553 //_______________________________________________________________________
554 Float_t AliEMCALRecoUtils::SmearClusterEnergy(const AliVCluster* cluster) 
555 {
556   //In case of MC analysis, smear energy to match resolution/calibration in real data
557   
558   if (!cluster) {
559     AliInfo("Cluster pointer null!");
560     return 0;
561   }
562   
563   Float_t energy    = cluster->E() ;
564   Float_t rdmEnergy = energy ;
565   if (fSmearClusterEnergy) {
566     rdmEnergy = fRandom.Gaus(energy,fSmearClusterParam[0] * TMath::Sqrt(energy) +
567                                     fSmearClusterParam[1] * energy +
568                                     fSmearClusterParam[2] );
569     AliDebug(2, Form("Energy: original %f, smeared %f\n", energy, rdmEnergy));
570   }
571   
572   return rdmEnergy;
573 }
574
575 //____________________________________________________________________________
576 Float_t AliEMCALRecoUtils::CorrectClusterEnergyLinearity(AliVCluster* cluster)
577 {
578   // Correct cluster energy from non linearity functions
579   
580   if (!cluster) {
581     AliInfo("Cluster pointer null!");
582     return 0;
583   }
584   
585   Float_t energy = cluster->E();
586
587   if (energy < 0.05) {
588     // Clusters with less than 50 MeV or negative are not possible
589     AliInfo(Form("Too Low Cluster energy!, E = %f < 0.05 GeV",energy));
590     return 0;
591   }
592   
593   switch (fNonLinearityFunction) 
594   {
595     case kPi0MC:
596     {
597       //Non-Linearity correction (from MC with function ([0]*exp(-[1]/E))+(([2]/([3]*2.*TMath::Pi())*exp(-(E-[4])^2/(2.*[3]^2)))))
598       //fNonLinearityParams[0] = 1.014;
599       //fNonLinearityParams[1] =-0.03329;
600       //fNonLinearityParams[2] =-0.3853;
601       //fNonLinearityParams[3] = 0.5423;
602       //fNonLinearityParams[4] =-0.4335;
603        energy *= (fNonLinearityParams[0]*exp(-fNonLinearityParams[1]/energy))+
604                   ((fNonLinearityParams[2]/(fNonLinearityParams[3]*2.*TMath::Pi())*
605                     exp(-(energy-fNonLinearityParams[4])*(energy-fNonLinearityParams[4])/(2.*fNonLinearityParams[3]*fNonLinearityParams[3]))));
606       break;
607     }
608      
609     case kPi0MCv2:
610     {
611       //Non-Linearity correction (from MC with function [0]/((x+[1])^[2]))+1;
612       //fNonLinearityParams[0] = 3.11111e-02;
613       //fNonLinearityParams[1] =-5.71666e-02; 
614       //fNonLinearityParams[2] = 5.67995e-01;      
615       
616       energy *= fNonLinearityParams[0]/TMath::Power(energy+fNonLinearityParams[1],fNonLinearityParams[2])+1;
617       break;
618     }
619     
620     case kPi0MCv3:
621     {
622       //Same as beam test corrected, change parameters
623       //fNonLinearityParams[0] =  9.81039e-01
624       //fNonLinearityParams[1] =  1.13508e-01;
625       //fNonLinearityParams[2] =  1.00173e+00; 
626       //fNonLinearityParams[3] =  9.67998e-02;
627       //fNonLinearityParams[4] =  2.19381e+02;
628       //fNonLinearityParams[5] =  6.31604e+01;
629       //fNonLinearityParams[6] =  1;
630       energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
631       
632       break;
633     }
634       
635       
636     case kPi0GammaGamma:
637     {
638       //Non-Linearity correction (from Olga Data with function p0+p1*exp(-p2*E))
639       //fNonLinearityParams[0] = 1.04;
640       //fNonLinearityParams[1] = -0.1445;
641       //fNonLinearityParams[2] = 1.046;
642       energy /= (fNonLinearityParams[0]+fNonLinearityParams[1]*exp(-fNonLinearityParams[2]*energy)); //Olga function
643       break;
644     }
645       
646     case kPi0GammaConversion:
647     {
648       //Non-Linearity correction (Nicolas from Dimitri Data with function C*[1-a*exp(-b*E)])
649       //fNonLinearityParams[0] = 0.139393/0.1349766;
650       //fNonLinearityParams[1] = 0.0566186;
651       //fNonLinearityParams[2] = 0.982133;
652       energy /= fNonLinearityParams[0]*(1-fNonLinearityParams[1]*exp(-fNonLinearityParams[2]*energy));
653       
654       break;
655     }
656       
657     case kBeamTest:
658     {
659       //From beam test, Alexei's results, for different ZS thresholds
660       //                        th=30 MeV; th = 45 MeV; th = 75 MeV
661       //fNonLinearityParams[0] = 1.007;      1.003;      1.002 
662       //fNonLinearityParams[1] = 0.894;      0.719;      0.797 
663       //fNonLinearityParams[2] = 0.246;      0.334;      0.358 
664       //Rescale the param[0] with 1.03
665       energy /= fNonLinearityParams[0]/(1+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]));
666       
667       break;
668     }
669       
670     case kBeamTestCorrected:
671     {
672       //From beam test, corrected for material between beam and EMCAL
673       //fNonLinearityParams[0] =  0.99078
674       //fNonLinearityParams[1] =  0.161499;
675       //fNonLinearityParams[2] =  0.655166; 
676       //fNonLinearityParams[3] =  0.134101;
677       //fNonLinearityParams[4] =  163.282;
678       //fNonLinearityParams[5] =  23.6904;
679       //fNonLinearityParams[6] =  0.978;
680         energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
681
682       break;
683     }
684      
685     case kBeamTestCorrectedv2:
686     {
687       //From beam test, corrected for material between beam and EMCAL
688       //fNonLinearityParams[0] =  0.983504;
689       //fNonLinearityParams[1] =  0.210106;
690       //fNonLinearityParams[2] =  0.897274;
691       //fNonLinearityParams[3] =  0.0829064;
692       //fNonLinearityParams[4] =  152.299;
693       //fNonLinearityParams[5] =  31.5028;
694       //fNonLinearityParams[6] =  0.968;
695       energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
696       
697       break;
698     }
699       
700     case kSDMv5:
701     {
702       //Based on fit to the MC/data using kNoCorrection on the data - utilizes symmetric decay method and kPi0MCv5(MC) - 28 Oct 2013
703       //fNonLinearityParams[0] =  1.0;
704       //fNonLinearityParams[1] =  6.64778e-02;
705       //fNonLinearityParams[2] =  1.570;
706       //fNonLinearityParams[3] =  9.67998e-02;
707       //fNonLinearityParams[4] =  2.19381e+02;
708       //fNonLinearityParams[5] =  6.31604e+01;
709       //fNonLinearityParams[6] =  1.01286;
710       energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5])))) * (0.964 + exp(-3.132-0.435*energy*2.0));
711       
712       break;
713     }
714       
715     case kPi0MCv5:
716     {
717       //Based on comparing MC truth information to the reconstructed energy of clusters. 
718       //fNonLinearityParams[0] =  1.0;
719       //fNonLinearityParams[1] =  6.64778e-02;
720       //fNonLinearityParams[2] =  1.570;
721       //fNonLinearityParams[3] =  9.67998e-02;
722       //fNonLinearityParams[4] =  2.19381e+02;
723       //fNonLinearityParams[5] =  6.31604e+01;
724       //fNonLinearityParams[6] =  1.01286;
725       energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
726       
727       break;
728     }
729       
730     case kNoCorrection:
731       AliDebug(2,"No correction on the energy\n");
732       break;
733       
734   }
735
736   return energy;
737 }
738
739 //__________________________________________________
740 void AliEMCALRecoUtils::InitNonLinearityParam()
741 {
742   //Initialising Non Linearity Parameters
743   
744   if (fNonLinearityFunction == kPi0MC) {
745     fNonLinearityParams[0] = 1.014;
746     fNonLinearityParams[1] = -0.03329;
747     fNonLinearityParams[2] = -0.3853;
748     fNonLinearityParams[3] = 0.5423;
749     fNonLinearityParams[4] = -0.4335;
750   }
751   
752   if (fNonLinearityFunction == kPi0MCv2) {
753     fNonLinearityParams[0] = 3.11111e-02;
754     fNonLinearityParams[1] =-5.71666e-02; 
755     fNonLinearityParams[2] = 5.67995e-01;      
756   }
757   
758   if (fNonLinearityFunction == kPi0MCv3) {
759     fNonLinearityParams[0] =  9.81039e-01;
760     fNonLinearityParams[1] =  1.13508e-01;
761     fNonLinearityParams[2] =  1.00173e+00; 
762     fNonLinearityParams[3] =  9.67998e-02;
763     fNonLinearityParams[4] =  2.19381e+02;
764     fNonLinearityParams[5] =  6.31604e+01;
765     fNonLinearityParams[6] =  1;
766   }
767   
768   if (fNonLinearityFunction == kPi0GammaGamma) {
769     fNonLinearityParams[0] = 1.04;
770     fNonLinearityParams[1] = -0.1445;
771     fNonLinearityParams[2] = 1.046;
772   }  
773
774   if (fNonLinearityFunction == kPi0GammaConversion) {
775     fNonLinearityParams[0] = 0.139393;
776     fNonLinearityParams[1] = 0.0566186;
777     fNonLinearityParams[2] = 0.982133;
778   }  
779
780   if (fNonLinearityFunction == kBeamTest) {
781     if (fNonLinearThreshold == 30) {
782       fNonLinearityParams[0] = 1.007; 
783       fNonLinearityParams[1] = 0.894; 
784       fNonLinearityParams[2] = 0.246; 
785     }
786     if (fNonLinearThreshold == 45) {
787       fNonLinearityParams[0] = 1.003; 
788       fNonLinearityParams[1] = 0.719; 
789       fNonLinearityParams[2] = 0.334; 
790     }
791     if (fNonLinearThreshold == 75) {
792       fNonLinearityParams[0] = 1.002; 
793       fNonLinearityParams[1] = 0.797; 
794       fNonLinearityParams[2] = 0.358; 
795     }
796   }
797
798   if (fNonLinearityFunction == kBeamTestCorrected) {
799     fNonLinearityParams[0] =  0.99078;
800     fNonLinearityParams[1] =  0.161499;
801     fNonLinearityParams[2] =  0.655166; 
802     fNonLinearityParams[3] =  0.134101;
803     fNonLinearityParams[4] =  163.282;
804     fNonLinearityParams[5] =  23.6904;
805     fNonLinearityParams[6] =  0.978;
806   }
807   
808   if (fNonLinearityFunction == kBeamTestCorrectedv2) {
809     fNonLinearityParams[0] =  0.983504;
810     fNonLinearityParams[1] =  0.210106;
811     fNonLinearityParams[2] =  0.897274;
812     fNonLinearityParams[3] =  0.0829064;
813     fNonLinearityParams[4] =  152.299;
814     fNonLinearityParams[5] =  31.5028;
815     fNonLinearityParams[6] =  0.968;
816   }
817
818   if (fNonLinearityFunction == kSDMv5) {
819     fNonLinearityParams[0] =  1.0;
820     fNonLinearityParams[1] =  6.64778e-02;
821     fNonLinearityParams[2] =  1.570;
822     fNonLinearityParams[3] =  9.67998e-02;
823     fNonLinearityParams[4] =  2.19381e+02;
824     fNonLinearityParams[5] =  6.31604e+01;
825     fNonLinearityParams[6] =  1.01286;
826   }
827
828   if (fNonLinearityFunction == kPi0MCv5) {
829     fNonLinearityParams[0] =  1.0;
830     fNonLinearityParams[1] =  6.64778e-02;
831     fNonLinearityParams[2] =  1.570;
832     fNonLinearityParams[3] =  9.67998e-02;
833     fNonLinearityParams[4] =  2.19381e+02;
834     fNonLinearityParams[5] =  6.31604e+01;
835     fNonLinearityParams[6] =  1.01286;
836   }
837
838 }
839
840 //_________________________________________________________
841 Float_t  AliEMCALRecoUtils::GetDepth(Float_t energy, 
842                                      Int_t iParticle, 
843                                      Int_t iSM) const 
844 {
845   //Calculate shower depth for a given cluster energy and particle type
846
847   // parameters 
848   Float_t x0    = 1.31;
849   Float_t ecr   = 8;
850   Float_t depth = 0;
851   Float_t arg   = energy*1000/ ecr; //Multiply energy by 1000 to transform to MeV
852   
853   switch ( iParticle )
854   {
855     case kPhoton:
856       if (arg < 1) 
857         depth = 0;
858       else
859         depth = x0 * (TMath::Log(arg) + 0.5); 
860       break;
861       
862     case kElectron:
863       if (arg < 1) 
864         depth = 0;
865       else
866         depth = x0 * (TMath::Log(arg) - 0.5); 
867       break;
868       
869     case kHadron:
870       // hadron 
871       // boxes anc. here
872       if (gGeoManager) {
873         gGeoManager->cd("ALIC_1/XEN1_1");
874         TGeoNode        *geoXEn1    = gGeoManager->GetCurrentNode();
875         TGeoNodeMatrix  *geoSM      = dynamic_cast<TGeoNodeMatrix *>(geoXEn1->GetDaughter(iSM));
876         if (geoSM) {
877           TGeoVolume      *geoSMVol   = geoSM->GetVolume(); 
878           TGeoShape       *geoSMShape = geoSMVol->GetShape();
879           TGeoBBox        *geoBox     = dynamic_cast<TGeoBBox *>(geoSMShape);
880           if (geoBox) depth = 0.5 * geoBox->GetDX()*2 ;
881           else AliFatal("Null GEANT box");
882         }
883         else AliFatal("NULL  GEANT node matrix");
884       }
885       else
886       {//electron
887         if (arg < 1) 
888           depth = 0;
889         else
890           depth = x0 * (TMath::Log(arg) - 0.5); 
891       }
892         
893       break;
894       
895     default://photon
896       if (arg < 1) 
897         depth = 0;
898       else
899         depth = x0 * (TMath::Log(arg) + 0.5);
900   }  
901   
902   return depth;
903 }
904
905 //____________________________________________________________________
906 void AliEMCALRecoUtils::GetMaxEnergyCell(const AliEMCALGeometry *geom, 
907                                          AliVCaloCells* cells, 
908                                          const AliVCluster* clu, 
909                                          Int_t  & absId,  
910                                          Int_t  & iSupMod, 
911                                          Int_t  & ieta, 
912                                          Int_t  & iphi, 
913                                          Bool_t & shared)
914 {
915   //For a given CaloCluster gets the absId of the cell 
916   //with maximum energy deposit.
917   
918   Double_t eMax        = -1.;
919   Double_t eCell       = -1.;
920   Float_t  fraction    = 1.;
921   Float_t  recalFactor = 1.;
922   Int_t    cellAbsId   = -1 ;
923
924   Int_t iTower  = -1;
925   Int_t iIphi   = -1;
926   Int_t iIeta   = -1;
927   Int_t iSupMod0= -1;
928
929   if (!clu) {
930     AliInfo("Cluster pointer null!");
931     absId=-1; iSupMod0=-1, ieta = -1; iphi = -1; shared = -1;
932     return;
933   }
934   
935   for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++) {
936     cellAbsId = clu->GetCellAbsId(iDig);
937     fraction  = clu->GetCellAmplitudeFraction(iDig);
938     //printf("a Cell %d, id, %d, amp %f, fraction %f\n",iDig,cellAbsId,cells->GetCellAmplitude(cellAbsId),fraction);
939     if (fraction < 1e-4) fraction = 1.; // in case unfolding is off
940     geom->GetCellIndex(cellAbsId,iSupMod,iTower,iIphi,iIeta); 
941     geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);
942     if (iDig==0) {
943       iSupMod0=iSupMod;
944     } else if (iSupMod0!=iSupMod)  {
945       shared = kTRUE;
946       //printf("AliEMCALRecoUtils::GetMaxEnergyCell() - SHARED CLUSTER\n");
947     }
948     if (!fCellsRecalibrated && IsRecalibrationOn()) {
949       recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
950     }
951     eCell  = cells->GetCellAmplitude(cellAbsId)*fraction*recalFactor;
952     //printf("b Cell %d, id, %d, amp %f, fraction %f\n",iDig,cellAbsId,eCell,fraction);
953     if (eCell > eMax) { 
954       eMax  = eCell; 
955       absId = cellAbsId;
956       //printf("\t new max: cell %d, e %f, ecell %f\n",maxId, eMax,eCell);
957     }
958   }// cell loop
959   
960   //Get from the absid the supermodule, tower and eta/phi numbers
961   geom->GetCellIndex(absId,iSupMod,iTower,iIphi,iIeta); 
962   //Gives SuperModule and Tower numbers
963   geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,
964                                          iIphi, iIeta,iphi,ieta); 
965   //printf("Max id %d, iSM %d, col %d, row %d\n",absId,iSupMod,ieta,iphi);
966   //printf("Max end---\n");
967 }
968
969 //______________________________________
970 void AliEMCALRecoUtils::InitParameters()
971 {
972   // Initialize data members with default values
973   
974   fParticleType = kPhoton;
975   fPosAlgo      = kUnchanged;
976   fW0           = 4.5;
977   
978   fNonLinearityFunction = kNoCorrection;
979   fNonLinearThreshold   = 30;
980   
981   fExoticCellFraction     = 0.97;
982   fExoticCellDiffTime     = 1e6;
983   fExoticCellMinAmplitude = 0.5;
984   
985   fAODFilterMask    = 128;
986   fAODHybridTracks  = kFALSE;
987   fAODTPCOnlyTracks = kTRUE;
988   
989   fCutEtaPhiSum      = kTRUE;
990   fCutEtaPhiSeparate = kFALSE;
991   
992   fCutR   = 0.05; 
993   fCutEta = 0.025; 
994   fCutPhi = 0.05;
995   
996   fClusterWindow = 100;
997   fMass          = 0.139;
998   
999   fStepSurface   = 20.;                      
1000   fStepCluster   = 5.;
1001   fTrackCutsType = kLooseCut;
1002   
1003   fCutMinTrackPt     = 0;
1004   fCutMinNClusterTPC = -1;
1005   fCutMinNClusterITS = -1;
1006   
1007   fCutMaxChi2PerClusterTPC  = 1e10;
1008   fCutMaxChi2PerClusterITS  = 1e10;
1009   
1010   fCutRequireTPCRefit     = kFALSE;
1011   fCutRequireITSRefit     = kFALSE;
1012   fCutAcceptKinkDaughters = kFALSE;
1013   
1014   fCutMaxDCAToVertexXY = 1e10;             
1015   fCutMaxDCAToVertexZ  = 1e10;              
1016   fCutDCAToVertex2D    = kFALSE;
1017   
1018   fCutRequireITSStandAlone = kFALSE; //MARCEL
1019   fCutRequireITSpureSA     = kFALSE; //Marcel
1020   
1021   //Misalignment matrices
1022   for (Int_t i = 0; i < 15 ; i++) 
1023   {
1024     fMisalTransShift[i] = 0.; 
1025     fMisalRotShift[i]   = 0.; 
1026   }
1027   
1028   //Non linearity
1029   for (Int_t i = 0; i < 7  ; i++) fNonLinearityParams[i] = 0.; 
1030   
1031   //For kBeamTestCorrectedv2 case, but default is no correction
1032   fNonLinearityParams[0] =  0.983504;
1033   fNonLinearityParams[1] =  0.210106;
1034   fNonLinearityParams[2] =  0.897274;
1035   fNonLinearityParams[3] =  0.0829064;
1036   fNonLinearityParams[4] =  152.299;
1037   fNonLinearityParams[5] =  31.5028;
1038   fNonLinearityParams[6] =  0.968;
1039   
1040   //Cluster energy smearing
1041   fSmearClusterEnergy   = kFALSE;
1042   fSmearClusterParam[0] = 0.07; // * sqrt E term
1043   fSmearClusterParam[1] = 0.00; // * E term
1044   fSmearClusterParam[2] = 0.00; // constant
1045 }
1046
1047 //_____________________________________________________
1048 void AliEMCALRecoUtils::InitEMCALRecalibrationFactors()
1049 {
1050   //Init EMCAL recalibration factors
1051   AliDebug(2,"AliCalorimeterUtils::InitEMCALRecalibrationFactors()");
1052   //In order to avoid rewriting the same histograms
1053   Bool_t oldStatus = TH1::AddDirectoryStatus();
1054   TH1::AddDirectory(kFALSE);
1055   
1056   fEMCALRecalibrationFactors = new TObjArray(12);
1057   for (int i = 0; i < 12; i++) 
1058     fEMCALRecalibrationFactors->Add(new TH2F(Form("EMCALRecalFactors_SM%d",i),
1059                                              Form("EMCALRecalFactors_SM%d",i),  48, 0, 48, 24, 0, 24));
1060   //Init the histograms with 1
1061   for (Int_t sm = 0; sm < 12; sm++) 
1062   {
1063     for (Int_t i = 0; i < 48; i++) 
1064     {
1065       for (Int_t j = 0; j < 24; j++) 
1066       {
1067         SetEMCALChannelRecalibrationFactor(sm,i,j,1.);
1068       }
1069     }
1070   }
1071   
1072   fEMCALRecalibrationFactors->SetOwner(kTRUE);
1073   fEMCALRecalibrationFactors->Compress();
1074   
1075   //In order to avoid rewriting the same histograms
1076   TH1::AddDirectory(oldStatus);    
1077 }
1078
1079 //_________________________________________________________
1080 void AliEMCALRecoUtils::InitEMCALTimeRecalibrationFactors()
1081 {
1082   //Init EMCAL recalibration factors
1083   AliDebug(2,"AliCalorimeterUtils::InitEMCALRecalibrationFactors()");
1084   //In order to avoid rewriting the same histograms
1085   Bool_t oldStatus = TH1::AddDirectoryStatus();
1086   TH1::AddDirectory(kFALSE);
1087   
1088   fEMCALTimeRecalibrationFactors = new TObjArray(4);
1089   for (int i = 0; i < 4; i++) 
1090     fEMCALTimeRecalibrationFactors->Add(new TH1F(Form("hAllTimeAvBC%d",i),
1091                                                  Form("hAllTimeAvBC%d",i),  
1092                                                  48*24*12,0.,48*24*12)          );
1093   //Init the histograms with 1
1094   for (Int_t bc = 0; bc < 4; bc++) 
1095   {
1096     for (Int_t i = 0; i < 48*24*12; i++) 
1097       SetEMCALChannelTimeRecalibrationFactor(bc,i,0.);
1098   }
1099   
1100   fEMCALTimeRecalibrationFactors->SetOwner(kTRUE);
1101   fEMCALTimeRecalibrationFactors->Compress();
1102   
1103   //In order to avoid rewriting the same histograms
1104   TH1::AddDirectory(oldStatus);    
1105 }
1106
1107 //____________________________________________________
1108 void AliEMCALRecoUtils::InitEMCALBadChannelStatusMap()
1109 {
1110   //Init EMCAL bad channels map
1111   AliDebug(2,"AliEMCALRecoUtils::InitEMCALBadChannelStatusMap()");
1112   //In order to avoid rewriting the same histograms
1113   Bool_t oldStatus = TH1::AddDirectoryStatus();
1114   TH1::AddDirectory(kFALSE);
1115   
1116   fEMCALBadChannelMap = new TObjArray(12);
1117   //TH2F * hTemp = new  TH2I("EMCALBadChannelMap","EMCAL SuperModule bad channel map", 48, 0, 48, 24, 0, 24);
1118   for (int i = 0; i < 12; i++) 
1119   {
1120     fEMCALBadChannelMap->Add(new TH2I(Form("EMCALBadChannelMap_Mod%d",i),Form("EMCALBadChannelMap_Mod%d",i), 48, 0, 48, 24, 0, 24));
1121   }
1122   
1123   fEMCALBadChannelMap->SetOwner(kTRUE);
1124   fEMCALBadChannelMap->Compress();
1125   
1126   //In order to avoid rewriting the same histograms
1127   TH1::AddDirectory(oldStatus);    
1128 }
1129
1130 //____________________________________________________________________________
1131 void AliEMCALRecoUtils::RecalibrateClusterEnergy(const AliEMCALGeometry* geom, 
1132                                                  AliVCluster * cluster, 
1133                                                  AliVCaloCells * cells, 
1134                                                  Int_t bc)
1135 {
1136   // Recalibrate the cluster energy and Time, considering the recalibration map 
1137   // and the energy of the cells and time that compose the cluster.
1138   // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1139   
1140   if (!cluster) {
1141     AliInfo("Cluster pointer null!");
1142     return;
1143   }  
1144   
1145   //Get the cluster number of cells and list of absId, check what kind of cluster do we have.
1146   UShort_t * index    = cluster->GetCellsAbsId() ;
1147   Double_t * fraction = cluster->GetCellsAmplitudeFraction() ;
1148   Int_t ncells = cluster->GetNCells();
1149   
1150   //Initialize some used variables
1151   Float_t energy = 0;
1152   Int_t   absId  =-1;
1153   Int_t   icol   =-1, irow =-1, imod=1;
1154   Float_t factor = 1, frac = 0;
1155   Int_t   absIdMax = -1;
1156   Float_t emax     = 0;
1157   
1158   //Loop on the cells, get the cell amplitude and recalibration factor, multiply and and to the new energy
1159   for (Int_t icell = 0; icell < ncells; icell++)
1160   {
1161     absId = index[icell];
1162     frac =  fraction[icell];
1163     if (frac < 1e-5) frac = 1; //in case of EMCAL, this is set as 0 since unfolding is off
1164     
1165     if (!fCellsRecalibrated && IsRecalibrationOn()) {
1166       // Energy  
1167       Int_t iTower = -1, iIphi = -1, iIeta = -1; 
1168       geom->GetCellIndex(absId,imod,iTower,iIphi,iIeta); 
1169       if (fEMCALRecalibrationFactors->GetEntries() <= imod) 
1170         continue;
1171       geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,irow,icol);      
1172       factor = GetEMCALChannelRecalibrationFactor(imod,icol,irow);
1173       
1174       AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - recalibrate cell: module %d, col %d, row %d, cell fraction %f,recalibration factor %f, cell energy %f\n",
1175                       imod,icol,irow,frac,factor,cells->GetCellAmplitude(absId)));
1176       
1177     } 
1178     
1179     energy += cells->GetCellAmplitude(absId)*factor*frac;
1180     
1181     if (emax < cells->GetCellAmplitude(absId)*factor*frac) {
1182       emax     = cells->GetCellAmplitude(absId)*factor*frac;
1183       absIdMax = absId;
1184     }
1185   }
1186   
1187   AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - Energy before %f, after %f \n",cluster->E(),energy));
1188
1189   cluster->SetE(energy);
1190
1191   // Recalculate time of cluster
1192   Double_t timeorg = cluster->GetTOF();
1193
1194   Double_t time = cells->GetCellTime(absIdMax);
1195   if (!fCellsRecalibrated && IsTimeRecalibrationOn())
1196     RecalibrateCellTime(absIdMax,bc,time);
1197
1198   cluster->SetTOF(time);
1199
1200   AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - Time before %f, after %f \n",timeorg,cluster->GetTOF()));
1201 }
1202
1203 //_____________________________________________________________
1204 void AliEMCALRecoUtils::RecalibrateCells(AliVCaloCells * cells,
1205                                          Int_t bc)
1206 {
1207   // Recalibrate the cells time and energy, considering the recalibration map and the energy 
1208   // of the cells that compose the cluster.
1209   // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1210
1211   if (!IsRecalibrationOn() && !IsTimeRecalibrationOn() && !IsBadChannelsRemovalSwitchedOn()) 
1212     return;
1213   
1214   if (!cells) {
1215     AliInfo("Cells pointer null!");
1216     return;
1217   }  
1218   
1219   Short_t  absId  =-1;
1220   Bool_t   accept = kFALSE;
1221   Float_t  ecell  = 0;
1222   Double_t tcell  = 0;
1223   Double_t ecellin = 0;
1224   Double_t tcellin = 0;
1225   Int_t  mclabel = -1;
1226   Double_t efrac = 0;
1227   
1228   Int_t nEMcell  = cells->GetNumberOfCells() ;  
1229   for (Int_t iCell = 0; iCell < nEMcell; iCell++) 
1230   { 
1231     cells->GetCell( iCell, absId, ecellin, tcellin, mclabel, efrac );
1232     
1233     accept = AcceptCalibrateCell(absId, bc, ecell ,tcell ,cells); 
1234     if (!accept) {
1235       ecell = 0;
1236       tcell = -1;
1237     }
1238     
1239     //Set new values
1240     cells->SetCell(iCell,absId,ecell, tcell, mclabel, efrac);
1241   }
1242
1243   fCellsRecalibrated = kTRUE;
1244 }
1245
1246 //_______________________________________________________________________________________________________
1247 void AliEMCALRecoUtils::RecalibrateCellTime(Int_t absId, Int_t bc, Double_t & celltime) const
1248 {
1249   // Recalibrate time of cell with absID  considering the recalibration map 
1250   // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1251   
1252   if (!fCellsRecalibrated && IsTimeRecalibrationOn() && bc >= 0) {
1253     celltime -= GetEMCALChannelTimeRecalibrationFactor(bc%4,absId)*1.e-9;    ;  
1254   }
1255 }
1256   
1257 //______________________________________________________________________________
1258 void AliEMCALRecoUtils::RecalculateClusterPosition(const AliEMCALGeometry *geom, 
1259                                                    AliVCaloCells* cells, 
1260                                                    AliVCluster* clu)
1261 {
1262   //For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1263   
1264   if (!clu) {
1265     AliInfo("Cluster pointer null!");
1266     return;
1267   }
1268     
1269   if      (fPosAlgo==kPosTowerGlobal) RecalculateClusterPositionFromTowerGlobal( geom, cells, clu);
1270   else if (fPosAlgo==kPosTowerIndex)  RecalculateClusterPositionFromTowerIndex ( geom, cells, clu);
1271   else    AliDebug(2,"Algorithm to recalculate position not selected, do nothing.");
1272 }  
1273
1274 //_____________________________________________________________________________________________
1275 void AliEMCALRecoUtils::RecalculateClusterPositionFromTowerGlobal(const AliEMCALGeometry *geom, 
1276                                                                   AliVCaloCells* cells, 
1277                                                                   AliVCluster* clu)
1278 {
1279   // For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1280   // The algorithm is a copy of what is done in AliEMCALRecPoint
1281   
1282   Double_t eCell       = 0.;
1283   Float_t  fraction    = 1.;
1284   Float_t  recalFactor = 1.;
1285   
1286   Int_t    absId   = -1;
1287   Int_t    iTower  = -1, iIphi  = -1, iIeta  = -1;
1288   Int_t    iSupModMax = -1, iSM=-1, iphi   = -1, ieta   = -1;
1289   Float_t  weight = 0.,  totalWeight=0.;
1290   Float_t  newPos[3] = {0,0,0};
1291   Double_t pLocal[3], pGlobal[3];
1292   Bool_t shared = kFALSE;
1293
1294   Float_t  clEnergy = clu->E(); //Energy already recalibrated previously
1295   if (clEnergy <= 0)
1296     return;
1297   GetMaxEnergyCell(geom, cells, clu, absId,  iSupModMax, ieta, iphi,shared);
1298   Double_t depth = GetDepth(clEnergy,fParticleType,iSupModMax) ;
1299   
1300   //printf("** Cluster energy %f, ncells %d, depth %f\n",clEnergy,clu->GetNCells(),depth);
1301   
1302   for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++) 
1303   {
1304     absId = clu->GetCellAbsId(iDig);
1305     fraction  = clu->GetCellAmplitudeFraction(iDig);
1306     if (fraction < 1e-4) fraction = 1.; // in case unfolding is off
1307     
1308     if (!fCellsRecalibrated) {
1309       geom->GetCellIndex(absId,iSM,iTower,iIphi,iIeta); 
1310       geom->GetCellPhiEtaIndexInSModule(iSM,iTower,iIphi, iIeta,iphi,ieta);      
1311       if (IsRecalibrationOn()) {
1312         recalFactor = GetEMCALChannelRecalibrationFactor(iSM,ieta,iphi);
1313       }
1314     }
1315     
1316     eCell  = cells->GetCellAmplitude(absId)*fraction*recalFactor;
1317     
1318     weight = GetCellWeight(eCell,clEnergy);
1319     totalWeight += weight;
1320     
1321     geom->RelPosCellInSModule(absId,depth,pLocal[0],pLocal[1],pLocal[2]);
1322     //printf("pLocal (%f,%f,%f), SM %d, absId %d\n",pLocal[0],pLocal[1],pLocal[2],iSupModMax,absId);
1323     geom->GetGlobal(pLocal,pGlobal,iSupModMax);
1324     //printf("pLocal (%f,%f,%f)\n",pGlobal[0],pGlobal[1],pGlobal[2]);
1325
1326     for (int i=0; i<3; i++ ) newPos[i] += (weight*pGlobal[i]);
1327   }// cell loop
1328   
1329   if (totalWeight>0) {
1330     for (int i=0; i<3; i++ )    newPos[i] /= totalWeight;
1331   }
1332     
1333   //Float_t pos[]={0,0,0};
1334   //clu->GetPosition(pos);
1335   //printf("OldPos  : %2.3f,%2.3f,%2.3f\n",pos[0],pos[1],pos[2]);
1336   //printf("NewPos  : %2.3f,%2.3f,%2.3f\n",newPos[0],newPos[1],newPos[2]);
1337   
1338   if (iSupModMax > 1) { //sector 1
1339     newPos[0] +=fMisalTransShift[3];//-=3.093; 
1340     newPos[1] +=fMisalTransShift[4];//+=6.82;
1341     newPos[2] +=fMisalTransShift[5];//+=1.635;
1342     //printf("   +    : %2.3f,%2.3f,%2.3f\n",fMisalTransShift[3],fMisalTransShift[4],fMisalTransShift[5]);
1343   } else { //sector 0
1344     newPos[0] +=fMisalTransShift[0];//+=1.134;
1345     newPos[1] +=fMisalTransShift[1];//+=8.2;
1346     newPos[2] +=fMisalTransShift[2];//+=1.197;
1347     //printf("   +    : %2.3f,%2.3f,%2.3f\n",fMisalTransShift[0],fMisalTransShift[1],fMisalTransShift[2]);
1348   }
1349   //printf("NewPos : %2.3f,%2.3f,%2.3f\n",newPos[0],newPos[1],newPos[2]);
1350
1351   clu->SetPosition(newPos);
1352 }  
1353
1354 //____________________________________________________________________________________________
1355 void AliEMCALRecoUtils::RecalculateClusterPositionFromTowerIndex(const AliEMCALGeometry *geom, 
1356                                                                  AliVCaloCells* cells, 
1357                                                                  AliVCluster* clu)
1358 {
1359   // For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1360   // The algorithm works with the tower indeces, averages the indeces and from them it calculates the global position
1361   
1362   Double_t eCell       = 1.;
1363   Float_t  fraction    = 1.;
1364   Float_t  recalFactor = 1.;
1365   
1366   Int_t absId   = -1;
1367   Int_t iTower  = -1;
1368   Int_t iIphi   = -1, iIeta   = -1;
1369   Int_t iSupMod = -1, iSupModMax = -1;
1370   Int_t iphi = -1, ieta =-1;
1371   Bool_t shared = kFALSE;
1372
1373   Float_t clEnergy = clu->E(); //Energy already recalibrated previously.
1374   
1375   if (clEnergy <= 0)
1376     return;
1377   GetMaxEnergyCell(geom, cells, clu, absId,  iSupModMax, ieta, iphi,shared);
1378   Float_t  depth = GetDepth(clEnergy,fParticleType,iSupMod) ;
1379
1380   Float_t weight = 0., weightedCol = 0., weightedRow = 0., totalWeight=0.;
1381   Bool_t areInSameSM = kTRUE; //exclude clusters with cells in different SMs for now
1382   Int_t startingSM = -1;
1383   
1384   for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++) 
1385   {
1386     absId = clu->GetCellAbsId(iDig);
1387     fraction  = clu->GetCellAmplitudeFraction(iDig);
1388     if (fraction < 1e-4) fraction = 1.; // in case unfolding is off
1389
1390     if      (iDig==0)  startingSM = iSupMod;
1391     else if (iSupMod != startingSM) areInSameSM = kFALSE;
1392
1393     eCell  = cells->GetCellAmplitude(absId);
1394     
1395     geom->GetCellIndex(absId,iSupMod,iTower,iIphi,iIeta); 
1396     geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);    
1397     
1398     if (!fCellsRecalibrated)
1399     {
1400       if (IsRecalibrationOn()) {
1401         recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1402       }
1403     }
1404     
1405     eCell  = cells->GetCellAmplitude(absId)*fraction*recalFactor;
1406     
1407     weight = GetCellWeight(eCell,clEnergy);
1408     if (weight < 0) weight = 0;
1409     totalWeight += weight;
1410     weightedCol += ieta*weight;
1411     weightedRow += iphi*weight;
1412     
1413     //printf("Max cell? cell %d, amplitude org %f, fraction %f, recalibration %f, amplitude new %f \n",cellAbsId, cells->GetCellAmplitude(cellAbsId), fraction, recalFactor, eCell) ;
1414   }// cell loop
1415     
1416   Float_t xyzNew[]={0.,0.,0.};
1417   if (areInSameSM == kTRUE) {
1418     //printf("In Same SM\n");
1419     weightedCol = weightedCol/totalWeight;
1420     weightedRow = weightedRow/totalWeight;
1421     geom->RecalculateTowerPosition(weightedRow, weightedCol, iSupModMax, depth, fMisalTransShift, fMisalRotShift, xyzNew); 
1422   } 
1423   else 
1424   {
1425     //printf("In Different SM\n");
1426     geom->RecalculateTowerPosition(iphi,        ieta,        iSupModMax, depth, fMisalTransShift, fMisalRotShift, xyzNew); 
1427   }
1428   
1429   clu->SetPosition(xyzNew);
1430 }
1431
1432 //___________________________________________________________________________________________
1433 void AliEMCALRecoUtils::RecalculateClusterDistanceToBadChannel(const AliEMCALGeometry * geom, 
1434                                                                AliVCaloCells* cells, 
1435                                                                AliVCluster * cluster)
1436 {           
1437   //re-evaluate distance to bad channel with updated bad map
1438   
1439   if (!fRecalDistToBadChannels) return;
1440   
1441   if (!cluster)
1442   {
1443     AliInfo("Cluster pointer null!");
1444     return;
1445   }  
1446   
1447   //Get channels map of the supermodule where the cluster is.
1448   Int_t absIdMax  = -1, iSupMod =-1, icolM = -1, irowM = -1;
1449   Bool_t shared = kFALSE;
1450   GetMaxEnergyCell(geom, cells, cluster, absIdMax,  iSupMod, icolM, irowM, shared);
1451   TH2D* hMap  = (TH2D*)fEMCALBadChannelMap->At(iSupMod);
1452
1453   Int_t dRrow, dRcol;  
1454   Float_t  minDist = 10000.;
1455   Float_t  dist    = 0.;
1456   
1457   //Loop on tower status map 
1458   for (Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++)
1459   {
1460     for (Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++)
1461     {
1462       //Check if tower is bad.
1463       if (hMap->GetBinContent(icol,irow)==0) continue;
1464       //printf("AliEMCALRecoUtils::RecalculateDistanceToBadChannels() - \n \t Bad channel in SM %d, col %d, row %d, \n \t Cluster max in col %d, row %d\n",
1465       //       iSupMod,icol, irow, icolM,irowM);
1466       
1467       dRrow=TMath::Abs(irowM-irow);
1468       dRcol=TMath::Abs(icolM-icol);
1469       dist=TMath::Sqrt(dRrow*dRrow+dRcol*dRcol);
1470       if (dist < minDist)
1471       {
1472         //printf("MIN DISTANCE TO BAD %2.2f\n",dist);
1473         minDist = dist;
1474       }
1475     }
1476   }
1477   
1478   //In case the cluster is shared by 2 SuperModules, need to check the map of the second Super Module
1479   if (shared) 
1480   {
1481     TH2D* hMap2 = 0;
1482     Int_t iSupMod2 = -1;
1483     
1484     //The only possible combinations are (0,1), (2,3) ... (8,9)
1485     if (iSupMod%2) iSupMod2 = iSupMod-1;
1486     else           iSupMod2 = iSupMod+1;
1487     hMap2  = (TH2D*)fEMCALBadChannelMap->At(iSupMod2);
1488     
1489     //Loop on tower status map of second super module
1490     for (Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++)
1491     {
1492       for (Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++)
1493       {
1494         //Check if tower is bad.
1495         if (hMap2->GetBinContent(icol,irow)==0) 
1496           continue;
1497         //printf("AliEMCALRecoUtils::RecalculateDistanceToBadChannels(shared) - \n \t Bad channel in SM %d, col %d, row %d \n \t Cluster max in SM %d, col %d, row %d\n",
1498         //     iSupMod2,icol, irow,iSupMod,icolM,irowM);
1499         dRrow=TMath::Abs(irow-irowM);
1500         
1501         if (iSupMod%2) {
1502           dRcol=TMath::Abs(icol-(AliEMCALGeoParams::fgkEMCALCols+icolM));
1503         } else {
1504           dRcol=TMath::Abs(AliEMCALGeoParams::fgkEMCALCols+icol-icolM);
1505         }                    
1506         
1507         dist=TMath::Sqrt(dRrow*dRrow+dRcol*dRcol);
1508         if (dist < minDist) minDist = dist;        
1509       }
1510     }
1511   }// shared cluster in 2 SuperModules
1512   
1513   AliDebug(2,Form("Max cluster cell (SM,col,row)=(%d %d %d) - Distance to Bad Channel %2.2f",iSupMod, icolM, irowM, minDist));
1514   cluster->SetDistanceToBadChannel(minDist);
1515 }
1516
1517 //__________________________________________________________________
1518 void AliEMCALRecoUtils::RecalculateClusterPID(AliVCluster * cluster)
1519 {           
1520   //re-evaluate identification parameters with bayesian
1521   
1522   if (!cluster) {
1523     AliInfo("Cluster pointer null!");
1524     return;
1525   }
1526   
1527   if (cluster->GetM02() != 0)
1528     fPIDUtils->ComputePID(cluster->E(),cluster->GetM02());
1529   
1530   Float_t pidlist[AliPID::kSPECIESCN+1];
1531   for (Int_t i = 0; i < AliPID::kSPECIESCN+1; i++) pidlist[i] = fPIDUtils->GetPIDFinal(i);
1532         
1533   cluster->SetPID(pidlist);
1534 }
1535
1536 //___________________________________________________________________________________________________________________
1537 void AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(const AliEMCALGeometry * geom, 
1538                                                                 AliVCaloCells* cells, 
1539                                                                 AliVCluster * cluster,
1540                                                                 Float_t & l0,   Float_t & l1,   
1541                                                                 Float_t & disp, Float_t & dEta, Float_t & dPhi,
1542                                                                 Float_t & sEta, Float_t & sPhi, Float_t & sEtaPhi)
1543 {
1544   // Calculates new center of gravity in the local EMCAL-module coordinates 
1545   // and tranfers into global ALICE coordinates
1546   // Calculates Dispersion and main axis
1547   
1548   if (!cluster) {
1549     AliInfo("Cluster pointer null!");
1550     return;
1551   }
1552     
1553   Double_t eCell       = 0.;
1554   Float_t  fraction    = 1.;
1555   Float_t  recalFactor = 1.;
1556
1557   Int_t    iSupMod = -1;
1558   Int_t    iTower  = -1;
1559   Int_t    iIphi   = -1;
1560   Int_t    iIeta   = -1;
1561   Int_t    iphi    = -1;
1562   Int_t    ieta    = -1;
1563   Double_t etai    = -1.;
1564   Double_t phii    = -1.;
1565   
1566   Int_t    nstat   = 0 ;
1567   Float_t  wtot    = 0.;
1568   Double_t w       = 0.;
1569   Double_t etaMean = 0.;
1570   Double_t phiMean = 0.;
1571   
1572   //Loop on cells, calculate the cluster energy, in case a cut on cell energy is added
1573   // and to check if the cluster is between 2 SM in eta
1574   Int_t   iSM0   = -1;
1575   Bool_t  shared = kFALSE;
1576   Float_t energy = 0;
1577   
1578   for (Int_t iDigit=0; iDigit < cluster->GetNCells(); iDigit++)
1579   {
1580     //Get from the absid the supermodule, tower and eta/phi numbers
1581     geom->GetCellIndex(cluster->GetCellAbsId(iDigit),iSupMod,iTower,iIphi,iIeta);
1582     geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi,iIeta, iphi,ieta);
1583     
1584     //Check if there are cells of different SM
1585     if      (iDigit == 0   ) iSM0 = iSupMod;
1586     else if (iSupMod!= iSM0) shared = kTRUE;
1587     
1588     //Get the cell energy, if recalibration is on, apply factors
1589     fraction  = cluster->GetCellAmplitudeFraction(iDigit);
1590     if (fraction < 1e-4) fraction = 1.; // in case unfolding is off
1591     
1592     if (IsRecalibrationOn()) {
1593       recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1594     }
1595     
1596     eCell  = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1597     
1598     energy += eCell;
1599     
1600   }//cell loop
1601   
1602   //Loop on cells
1603   for (Int_t iDigit=0; iDigit < cluster->GetNCells(); iDigit++) 
1604   {
1605     //Get from the absid the supermodule, tower and eta/phi numbers
1606     geom->GetCellIndex(cluster->GetCellAbsId(iDigit),iSupMod,iTower,iIphi,iIeta); 
1607     geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi,iIeta, iphi,ieta);        
1608     
1609     //Get the cell energy, if recalibration is on, apply factors
1610     fraction  = cluster->GetCellAmplitudeFraction(iDigit);
1611     if (fraction < 1e-4) fraction = 1.; // in case unfolding is off
1612     
1613     if (!fCellsRecalibrated) {
1614       if (IsRecalibrationOn()) {
1615         recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1616       }
1617     }
1618     
1619     eCell  = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1620     
1621     // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
1622     // C Side impair SM, nSupMod%2=1; A side pair SM, nSupMod%2=0
1623     if (shared && iSupMod%2) ieta+=AliEMCALGeoParams::fgkEMCALCols;
1624     
1625     if (cluster->E() > 0 && eCell > 0) {
1626       w  = GetCellWeight(eCell,cluster->E());
1627       
1628       etai=(Double_t)ieta;
1629       phii=(Double_t)iphi;  
1630       
1631       if (w > 0.0) {
1632         wtot += w ;
1633         nstat++;            
1634         //Shower shape
1635         sEta     += w * etai * etai ;
1636         etaMean  += w * etai ;
1637         sPhi     += w * phii * phii ;
1638         phiMean  += w * phii ; 
1639         sEtaPhi  += w * etai * phii ; 
1640       }
1641     } else
1642       AliError(Form("Wrong energy %f and/or amplitude %f\n", eCell, cluster->E()));
1643   }//cell loop
1644   
1645   //Normalize to the weight  
1646   if (wtot > 0) {
1647     etaMean /= wtot ;
1648     phiMean /= wtot ;
1649   } else
1650     AliError(Form("Wrong weight %f\n", wtot));
1651   
1652   //Calculate dispersion  
1653   for (Int_t iDigit=0; iDigit < cluster->GetNCells(); iDigit++) 
1654   {
1655     //Get from the absid the supermodule, tower and eta/phi numbers
1656     geom->GetCellIndex(cluster->GetCellAbsId(iDigit),iSupMod,iTower,iIphi,iIeta); 
1657     geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi,iIeta, iphi,ieta);
1658     
1659     //Get the cell energy, if recalibration is on, apply factors
1660     fraction  = cluster->GetCellAmplitudeFraction(iDigit);
1661     if (fraction < 1e-4) fraction = 1.; // in case unfolding is off
1662     if (IsRecalibrationOn()) {
1663       recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1664     }
1665     eCell  = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1666     
1667     // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
1668     // C Side impair SM, nSupMod%2=1; A side pair SM, nSupMod%2=0
1669     if (shared && iSupMod%2) ieta+=AliEMCALGeoParams::fgkEMCALCols;
1670     
1671     if (cluster->E() > 0 && eCell > 0) {
1672       w  = GetCellWeight(eCell,cluster->E());
1673       
1674       etai=(Double_t)ieta;
1675       phii=(Double_t)iphi;    
1676       if (w > 0.0) { 
1677         disp +=  w *((etai-etaMean)*(etai-etaMean)+(phii-phiMean)*(phii-phiMean)); 
1678         dEta +=  w * (etai-etaMean)*(etai-etaMean) ; 
1679         dPhi +=  w * (phii-phiMean)*(phii-phiMean) ; 
1680       }
1681     } else
1682       AliError(Form("Wrong energy %f and/or amplitude %f\n", eCell, cluster->E()));
1683   }// cell loop
1684   
1685   //Normalize to the weigth and set shower shape parameters
1686   if (wtot > 0 && nstat > 1) {
1687     disp    /= wtot ;
1688     dEta    /= wtot ;
1689     dPhi    /= wtot ;
1690     sEta    /= wtot ;
1691     sPhi    /= wtot ;
1692     sEtaPhi /= wtot ;
1693     
1694     sEta    -= etaMean * etaMean ;
1695     sPhi    -= phiMean * phiMean ;
1696     sEtaPhi -= etaMean * phiMean ;
1697     
1698     l0 = (0.5 * (sEta + sPhi) + TMath::Sqrt( 0.25 * (sEta - sPhi) * (sEta - sPhi) + sEtaPhi * sEtaPhi ));
1699     l1 = (0.5 * (sEta + sPhi) - TMath::Sqrt( 0.25 * (sEta - sPhi) * (sEta - sPhi) + sEtaPhi * sEtaPhi ));
1700   } else {
1701     l0   = 0. ;
1702     l1   = 0. ;
1703     dEta = 0. ; dPhi = 0. ; disp    = 0. ;
1704     sEta = 0. ; sPhi = 0. ; sEtaPhi = 0. ;
1705   }  
1706 }
1707
1708 //____________________________________________________________________________________________
1709 void AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(const AliEMCALGeometry * geom, 
1710                                                                 AliVCaloCells* cells, 
1711                                                                 AliVCluster * cluster)
1712 {
1713   // Calculates new center of gravity in the local EMCAL-module coordinates 
1714   // and tranfers into global ALICE coordinates
1715   // Calculates Dispersion and main axis and puts them into the cluster
1716   
1717   Float_t l0   = 0., l1   = 0.;
1718   Float_t disp = 0., dEta = 0., dPhi    = 0.; 
1719   Float_t sEta = 0., sPhi = 0., sEtaPhi = 0.;
1720   
1721   AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(geom,cells,cluster,l0,l1,disp,
1722                                                              dEta, dPhi, sEta, sPhi, sEtaPhi);
1723   
1724   cluster->SetM02(l0);
1725   cluster->SetM20(l1);
1726   if (disp > 0. ) cluster->SetDispersion(TMath::Sqrt(disp)) ;
1727   
1728
1729
1730 //____________________________________________________________________________
1731 void AliEMCALRecoUtils::FindMatches(AliVEvent *event,
1732                                     TObjArray * clusterArr,  
1733                                     const AliEMCALGeometry *geom)
1734 {
1735   //This function should be called before the cluster loop
1736   //Before call this function, please recalculate the cluster positions
1737   //Given the input event, loop over all the tracks, select the closest cluster as matched with fCutR
1738   //Store matched cluster indexes and residuals
1739   
1740   fMatchedTrackIndex  ->Reset();
1741   fMatchedClusterIndex->Reset();
1742   fResidualPhi->Reset();
1743   fResidualEta->Reset();
1744   
1745   fMatchedTrackIndex  ->Set(1000);
1746   fMatchedClusterIndex->Set(1000);
1747   fResidualPhi->Set(1000);
1748   fResidualEta->Set(1000);
1749   
1750   AliESDEvent* esdevent = dynamic_cast<AliESDEvent*> (event);
1751   AliAODEvent* aodevent = dynamic_cast<AliAODEvent*> (event);
1752   
1753   // init the magnetic field if not already on
1754   if (!TGeoGlobalMagField::Instance()->GetField()) {
1755     if (!event->InitMagneticField())
1756     {
1757       AliInfo("Mag Field not initialized, null esd/aod evetn pointers");
1758     }
1759   } // Init mag field
1760   
1761   if (esdevent) {
1762     UInt_t mask1 = esdevent->GetESDRun()->GetDetectorsInDAQ();
1763     UInt_t mask2 = esdevent->GetESDRun()->GetDetectorsInReco();
1764     Bool_t desc1 = (mask1 >> 3) & 0x1;
1765     Bool_t desc2 = (mask2 >> 3) & 0x1;
1766     if (desc1==0 || desc2==0) { 
1767 //       AliError(Form("TPC not in DAQ/RECO: %u (%u)/%u (%u)", 
1768 //       mask1, esdevent->GetESDRun()->GetDetectorsInReco(),
1769 //       mask2, esdevent->GetESDRun()->GetDetectorsInDAQ()));
1770       fITSTrackSA=kTRUE;
1771     }
1772   }
1773
1774   TObjArray *clusterArray = 0x0;
1775   if (!clusterArr) {
1776     clusterArray = new TObjArray(event->GetNumberOfCaloClusters());
1777     for (Int_t icl=0; icl<event->GetNumberOfCaloClusters(); icl++) 
1778     {
1779       AliVCluster *cluster = (AliVCluster*) event->GetCaloCluster(icl);
1780       if (geom && !IsGoodCluster(cluster,geom,(AliVCaloCells*)event->GetEMCALCells())) 
1781         continue;
1782       clusterArray->AddAt(cluster,icl);
1783     }
1784   }
1785   
1786   Int_t    matched=0;
1787   Double_t cv[21];
1788   for (Int_t i=0; i<21;i++) cv[i]=0;
1789   for (Int_t itr=0; itr<event->GetNumberOfTracks(); itr++)
1790   {
1791     AliExternalTrackParam *trackParam = 0;
1792     
1793     //If the input event is ESD, the starting point for extrapolation is TPCOut, if available, or TPCInner 
1794     AliESDtrack *esdTrack = 0;
1795     AliAODTrack *aodTrack = 0;
1796     if (esdevent) {
1797       esdTrack = esdevent->GetTrack(itr);
1798       if (!esdTrack) continue;
1799       if (!IsAccepted(esdTrack)) continue;
1800       if (esdTrack->Pt()<fCutMinTrackPt) continue;
1801       Double_t phi = esdTrack->Phi()*TMath::RadToDeg();
1802       if (TMath::Abs(esdTrack->Eta())>0.9 || phi <= 10 || phi >= 250 ) continue;
1803       if (!fITSTrackSA)
1804         trackParam =  const_cast<AliExternalTrackParam*>(esdTrack->GetInnerParam());  // if TPC Available
1805       else
1806         trackParam =  new AliExternalTrackParam(*esdTrack); // If ITS Track Standing alone              
1807     }
1808     
1809     //If the input event is AOD, the starting point for extrapolation is at vertex
1810     //AOD tracks are selected according to its filterbit.
1811     else if (aodevent) {
1812       aodTrack = dynamic_cast<AliAODTrack*>(aodevent->GetTrack(itr));
1813       if(!aodTrack) AliFatal("Not a standard AOD");
1814       if (!aodTrack) continue;
1815             
1816       if (fAODTPCOnlyTracks) { // Match with TPC only tracks, default from May 2013, before filter bit 32
1817         //printf("Match with TPC only tracks, accept? %d, test bit 128 <%d> \n", aodTrack->IsTPCOnly(), aodTrack->TestFilterMask(128));
1818         if (!aodTrack->IsTPCConstrained()) continue ;
1819       } else if (fAODHybridTracks) { // Match with hybrid tracks
1820         //printf("Match with Hybrid tracks, accept? %d \n", aodTrack->IsHybridGlobalConstrainedGlobal());
1821         if (!aodTrack->IsHybridGlobalConstrainedGlobal()) continue ;
1822       } else { // Match with tracks on a mask
1823         //printf("Match with tracks having filter bit mask %d, accept? %d \n",fAODFilterMask,aodTrack->TestFilterMask(fAODFilterMask));
1824         if (!aodTrack->TestFilterMask(fAODFilterMask) ) continue; //Select AOD tracks
1825       }
1826       
1827       if (aodTrack->Pt()<fCutMinTrackPt) continue;
1828
1829       Double_t phi = aodTrack->Phi()*TMath::RadToDeg();
1830       if (TMath::Abs(aodTrack->Eta())>0.9 || phi <= 10 || phi >= 250 ) 
1831         continue;
1832       Double_t pos[3],mom[3];
1833       aodTrack->GetXYZ(pos);
1834       aodTrack->GetPxPyPz(mom);
1835       AliDebug(5,Form("aod track: i=%d | pos=(%5.4f,%5.4f,%5.4f) | mom=(%5.4f,%5.4f,%5.4f) | charge=%d\n",itr,pos[0],pos[1],pos[2],mom[0],mom[1],mom[2],aodTrack->Charge()));
1836
1837       trackParam= new AliExternalTrackParam(pos,mom,cv,aodTrack->Charge());
1838     }
1839     
1840     //Return if the input data is not "AOD" or "ESD"
1841     else {
1842       printf("Wrong input data type! Should be \"AOD\" or \"ESD\"\n");
1843       if (clusterArray) {
1844         clusterArray->Clear();
1845         delete clusterArray;
1846       }
1847       return;
1848     }
1849     
1850     if (!trackParam) continue;
1851
1852     //Extrapolate the track to EMCal surface
1853     AliExternalTrackParam emcalParam(*trackParam);
1854     Float_t eta, phi, pt;
1855     if (!ExtrapolateTrackToEMCalSurface(&emcalParam, fEMCalSurfaceDistance, fMass, fStepSurface, eta, phi, pt)) {
1856       if (aodevent && trackParam) delete trackParam;
1857       if (fITSTrackSA && trackParam) delete trackParam;
1858       continue;
1859     }
1860
1861     if (TMath::Abs(eta)>0.75 || (phi) < 70*TMath::DegToRad() || (phi) > 190*TMath::DegToRad()) {
1862       if (aodevent && trackParam) delete trackParam;
1863       if (fITSTrackSA && trackParam) delete trackParam;
1864       continue;
1865     }
1866
1867     //Find matched clusters
1868     Int_t index = -1;
1869     Float_t dEta = -999, dPhi = -999;
1870     if (!clusterArr) {
1871       index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArray, dEta, dPhi);  
1872     } else {
1873       index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArr, dEta, dPhi);  
1874     }  
1875     
1876     if (index>-1) {
1877       fMatchedTrackIndex   ->AddAt(itr,matched);
1878       fMatchedClusterIndex ->AddAt(index,matched);
1879       fResidualEta         ->AddAt(dEta,matched);
1880       fResidualPhi         ->AddAt(dPhi,matched);
1881       matched++;
1882     }
1883     if (aodevent && trackParam) delete trackParam;
1884     if (fITSTrackSA && trackParam) delete trackParam;
1885   }//track loop
1886
1887   if (clusterArray) {
1888     clusterArray->Clear();
1889     delete clusterArray;
1890   }
1891   
1892   AliDebug(2,Form("Number of matched pairs = %d !\n",matched));
1893   
1894   fMatchedTrackIndex   ->Set(matched);
1895   fMatchedClusterIndex ->Set(matched);
1896   fResidualPhi         ->Set(matched);
1897   fResidualEta         ->Set(matched);
1898 }
1899
1900 //________________________________________________________________________________
1901 Int_t AliEMCALRecoUtils::FindMatchedClusterInEvent(const AliESDtrack *track, 
1902                                                    const AliVEvent *event, 
1903                                                    const AliEMCALGeometry *geom, 
1904                                                    Float_t &dEta, Float_t &dPhi)
1905 {
1906   //
1907   // This function returns the index of matched cluster to input track
1908   // Returns -1 if no match is found
1909   Int_t index = -1;
1910   Double_t phiV = track->Phi()*TMath::RadToDeg();
1911   if (TMath::Abs(track->Eta())>0.9 || phiV <= 10 || phiV >= 250 ) return index;
1912   AliExternalTrackParam *trackParam = 0;
1913   if (!fITSTrackSA)
1914     trackParam = const_cast<AliExternalTrackParam*>(track->GetInnerParam());  // If TPC
1915   else
1916     trackParam = new AliExternalTrackParam(*track);
1917     
1918   if (!trackParam) return index;
1919   AliExternalTrackParam emcalParam(*trackParam);
1920   Float_t eta, phi, pt;
1921
1922   if (!ExtrapolateTrackToEMCalSurface(&emcalParam, fEMCalSurfaceDistance, fMass, fStepSurface, eta, phi, pt))   {
1923     if (fITSTrackSA) delete trackParam;
1924     return index;
1925   }
1926   if (TMath::Abs(eta)>0.75 || (phi) < 70*TMath::DegToRad() || (phi) > 190*TMath::DegToRad()) {
1927     if (fITSTrackSA) delete trackParam;
1928     return index;
1929   }
1930   
1931   TObjArray *clusterArr = new TObjArray(event->GetNumberOfCaloClusters());
1932
1933   for (Int_t icl=0; icl<event->GetNumberOfCaloClusters(); icl++)
1934   {
1935     AliVCluster *cluster = (AliVCluster*) event->GetCaloCluster(icl);
1936     if (geom && !IsGoodCluster(cluster,geom,(AliVCaloCells*)event->GetEMCALCells())) continue;
1937     clusterArr->AddAt(cluster,icl);
1938   }
1939
1940   index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArr, dEta, dPhi);  
1941   clusterArr->Clear();
1942   delete clusterArr;
1943   if (fITSTrackSA) delete trackParam;
1944
1945   return index;
1946 }
1947
1948 //_______________________________________________________________________________________________
1949 Int_t  AliEMCALRecoUtils::FindMatchedClusterInClusterArr(const AliExternalTrackParam *emcalParam, 
1950                                                          AliExternalTrackParam *trkParam, 
1951                                                          const TObjArray * clusterArr, 
1952                                                          Float_t &dEta, Float_t &dPhi)
1953 {
1954   // Find matched cluster in array
1955   
1956   dEta=-999, dPhi=-999;
1957   Float_t dRMax = fCutR, dEtaMax=fCutEta, dPhiMax=fCutPhi;
1958   Int_t index = -1;
1959   Float_t tmpEta=-999, tmpPhi=-999;
1960
1961   Double_t exPos[3] = {0.,0.,0.};
1962   if (!emcalParam->GetXYZ(exPos)) return index;
1963
1964   Float_t clsPos[3] = {0.,0.,0.};
1965   for (Int_t icl=0; icl<clusterArr->GetEntriesFast(); icl++)
1966   {
1967     AliVCluster *cluster = dynamic_cast<AliVCluster*> (clusterArr->At(icl)) ;
1968     if (!cluster || !cluster->IsEMCAL()) continue;
1969     cluster->GetPosition(clsPos);
1970     Double_t dR = TMath::Sqrt(TMath::Power(exPos[0]-clsPos[0],2)+TMath::Power(exPos[1]-clsPos[1],2)+TMath::Power(exPos[2]-clsPos[2],2));
1971     if (dR > fClusterWindow) continue;
1972     
1973     AliExternalTrackParam trkPamTmp (*trkParam);//Retrieve the starting point every time before the extrapolation
1974     if (!ExtrapolateTrackToCluster(&trkPamTmp, cluster, fMass, fStepCluster, tmpEta, tmpPhi)) continue;
1975     if (fCutEtaPhiSum) {
1976       Float_t tmpR=TMath::Sqrt(tmpEta*tmpEta + tmpPhi*tmpPhi);
1977       if (tmpR<dRMax) {
1978         dRMax=tmpR;
1979         dEtaMax=tmpEta;
1980         dPhiMax=tmpPhi;
1981         index=icl;
1982       }
1983     } else if (fCutEtaPhiSeparate) {
1984       if (TMath::Abs(tmpEta)<TMath::Abs(dEtaMax) && TMath::Abs(tmpPhi)<TMath::Abs(dPhiMax)) {
1985         dEtaMax = tmpEta;
1986         dPhiMax = tmpPhi;
1987         index=icl;
1988       }
1989     } else {
1990       printf("Error: please specify your cut criteria\n");
1991       printf("To cut on sqrt(dEta^2+dPhi^2), use: SwitchOnCutEtaPhiSum()\n");
1992       printf("To cut on dEta and dPhi separately, use: SwitchOnCutEtaPhiSeparate()\n");
1993       return index;
1994     }
1995   }
1996
1997   dEta=dEtaMax;
1998   dPhi=dPhiMax;
1999
2000   return index;
2001 }
2002
2003 //------------------------------------------------------------------------------------
2004 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToEMCalSurface(AliVTrack *track,
2005                                                          Double_t emcalR, Double_t mass,
2006                                                          Double_t step, Double_t minpt,
2007                                                          Bool_t useMassForTracking)
2008
2009   // Extrapolate track to EMCAL surface
2010
2011   track->SetTrackPhiEtaPtOnEMCal(-999, -999, -999);
2012
2013   if (track->Pt()<minpt)
2014     return kFALSE;
2015
2016   Double_t phi = track->Phi()*TMath::RadToDeg();
2017   if (TMath::Abs(track->Eta())>0.9 || phi <= 10 || phi >= 250) 
2018     return kFALSE;
2019
2020   AliESDtrack *esdt = dynamic_cast<AliESDtrack*>(track);
2021   AliAODTrack *aodt = 0;
2022   if (!esdt) {
2023     aodt = dynamic_cast<AliAODTrack*>(track);
2024     if (!aodt)
2025       return kFALSE;
2026   }
2027
2028   // Select the mass hypothesis
2029   if ( mass < 0 )
2030   {
2031     Bool_t onlyTPC = kFALSE;
2032     if ( mass == -99 ) onlyTPC=kTRUE;
2033     
2034     if (esdt)
2035     {
2036       if ( useMassForTracking ) mass = esdt->GetMassForTracking();
2037       else                      mass = esdt->GetMass(onlyTPC);
2038     }
2039     else
2040     {
2041       if ( useMassForTracking ) mass = aodt->GetMassForTracking();
2042       else                      mass = aodt->M();
2043     }
2044   }
2045
2046   AliExternalTrackParam *trackParam = 0;
2047   if (esdt) {
2048     const AliExternalTrackParam *in = esdt->GetInnerParam();
2049     if (!in)
2050       return kFALSE;
2051     trackParam = new AliExternalTrackParam(*in);
2052   } else {
2053     Double_t xyz[3] = {0}, pxpypz[3] = {0}, cv[21] = {0};
2054     aodt->PxPyPz(pxpypz);  
2055     aodt->XvYvZv(xyz);
2056     aodt->GetCovarianceXYZPxPyPz(cv);  
2057     trackParam = new AliExternalTrackParam(xyz,pxpypz,cv,aodt->Charge());
2058   }
2059   if (!trackParam)
2060     return kFALSE;
2061
2062   Float_t etaout=-999, phiout=-999, ptout=-999;
2063   Bool_t ret = ExtrapolateTrackToEMCalSurface(trackParam, 
2064                                               emcalR,
2065                                               mass,
2066                                               step,
2067                                               etaout, 
2068                                               phiout,
2069                                               ptout);
2070   delete trackParam;
2071   if (!ret)
2072     return kFALSE;
2073   if (TMath::Abs(etaout)>0.75 || (phiout<70*TMath::DegToRad()) || (phiout>190*TMath::DegToRad()))
2074     return kFALSE;
2075   track->SetTrackPhiEtaPtOnEMCal(phiout, etaout, ptout);
2076   return kTRUE;
2077 }
2078
2079
2080 //------------------------------------------------------------------------------------
2081 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToEMCalSurface(AliExternalTrackParam *trkParam, 
2082                                                          Double_t emcalR,
2083                                                          Double_t mass, 
2084                                                          Double_t step, 
2085                                                          Float_t &eta, 
2086                                                          Float_t &phi,
2087                                                          Float_t &pt)
2088 {
2089   //Extrapolate track to EMCAL surface
2090   
2091   eta = -999, phi = -999, pt = -999;
2092   if (!trkParam) return kFALSE;
2093   if (!AliTrackerBase::PropagateTrackToBxByBz(trkParam, emcalR, mass, step, kTRUE, 0.8, -1)) return kFALSE;
2094   Double_t trkPos[3] = {0.,0.,0.};
2095   if (!trkParam->GetXYZ(trkPos)) return kFALSE;
2096   TVector3 trkPosVec(trkPos[0],trkPos[1],trkPos[2]);
2097   eta = trkPosVec.Eta();
2098   phi = trkPosVec.Phi();
2099   pt = trkParam->Pt();
2100   if (phi<0)
2101     phi += 2*TMath::Pi();
2102
2103   return kTRUE;
2104 }
2105
2106 //-----------------------------------------------------------------------------------
2107 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToPosition(AliExternalTrackParam *trkParam, 
2108                                                      const Float_t *clsPos, 
2109                                                      Double_t mass, 
2110                                                      Double_t step, 
2111                                                      Float_t &tmpEta, 
2112                                                      Float_t &tmpPhi)
2113 {
2114   //
2115   //Return the residual by extrapolating a track param to a global position
2116   //
2117   tmpEta = -999;
2118   tmpPhi = -999;
2119   if (!trkParam) return kFALSE;
2120   Double_t trkPos[3] = {0.,0.,0.};
2121   TVector3 vec(clsPos[0],clsPos[1],clsPos[2]);
2122   Double_t alpha =  ((int)(vec.Phi()*TMath::RadToDeg()/20)+0.5)*20*TMath::DegToRad();
2123   vec.RotateZ(-alpha); //Rotate the cluster to the local extrapolation coordinate system
2124   if (!AliTrackerBase::PropagateTrackToBxByBz(trkParam, vec.X(), mass, step,kTRUE, 0.8, -1)) return kFALSE;
2125   if (!trkParam->GetXYZ(trkPos)) return kFALSE; //Get the extrapolated global position
2126
2127   TVector3 clsPosVec(clsPos[0],clsPos[1],clsPos[2]);
2128   TVector3 trkPosVec(trkPos[0],trkPos[1],trkPos[2]);
2129
2130   // track cluster matching
2131   tmpPhi = clsPosVec.DeltaPhi(trkPosVec);    // tmpPhi is between -pi and pi
2132   tmpEta = clsPosVec.Eta()-trkPosVec.Eta();
2133
2134   return kTRUE;
2135 }
2136
2137 //----------------------------------------------------------------------------------
2138 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToCluster(AliExternalTrackParam *trkParam, 
2139                                                     const AliVCluster *cluster, 
2140                                                     Double_t mass, 
2141                                                     Double_t step, 
2142                                                     Float_t &tmpEta, 
2143                                                     Float_t &tmpPhi)
2144 {
2145   //
2146   //Return the residual by extrapolating a track param to a cluster
2147   //
2148   tmpEta = -999;
2149   tmpPhi = -999;
2150   if (!cluster || !trkParam) 
2151     return kFALSE;
2152
2153   Float_t clsPos[3] = {0.,0.,0.};
2154   cluster->GetPosition(clsPos);
2155
2156   return ExtrapolateTrackToPosition(trkParam, clsPos, mass, step, tmpEta, tmpPhi);
2157 }
2158
2159 //---------------------------------------------------------------------------------
2160 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToCluster(AliExternalTrackParam *trkParam, 
2161                                                     const AliVCluster *cluster, 
2162                                                     Float_t &tmpEta, 
2163                                                     Float_t &tmpPhi)
2164 {
2165   //
2166   //Return the residual by extrapolating a track param to a clusterfStepCluster
2167   //
2168
2169   return ExtrapolateTrackToCluster(trkParam, cluster, fMass, fStepCluster, tmpEta, tmpPhi);
2170 }
2171
2172 //_______________________________________________________________________
2173 void AliEMCALRecoUtils::GetMatchedResiduals(Int_t clsIndex, 
2174                                             Float_t &dEta, Float_t &dPhi)
2175 {
2176   //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2177   //Get the residuals dEta and dPhi for this cluster to the closest track
2178   //Works with ESDs and AODs
2179
2180   if (FindMatchedPosForCluster(clsIndex) >= 999) {
2181     AliDebug(2,"No matched tracks found!\n");
2182     dEta=999.;
2183     dPhi=999.;
2184     return;
2185   }
2186   dEta = fResidualEta->At(FindMatchedPosForCluster(clsIndex));
2187   dPhi = fResidualPhi->At(FindMatchedPosForCluster(clsIndex));
2188 }
2189
2190 //______________________________________________________________________________________________
2191 void AliEMCALRecoUtils::GetMatchedClusterResiduals(Int_t trkIndex, Float_t &dEta, Float_t &dPhi)
2192 {
2193   //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2194   //Get the residuals dEta and dPhi for this track to the closest cluster
2195   //Works with ESDs and AODs
2196
2197   if (FindMatchedPosForTrack(trkIndex) >= 999) {
2198     AliDebug(2,"No matched cluster found!\n");
2199     dEta=999.;
2200     dPhi=999.;
2201     return;
2202   }
2203   dEta = fResidualEta->At(FindMatchedPosForTrack(trkIndex));
2204   dPhi = fResidualPhi->At(FindMatchedPosForTrack(trkIndex));
2205 }
2206
2207 //__________________________________________________________
2208 Int_t AliEMCALRecoUtils::GetMatchedTrackIndex(Int_t clsIndex)
2209 {
2210   //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2211   //Get the index of matched track to this cluster
2212   //Works with ESDs and AODs
2213   
2214   if (IsClusterMatched(clsIndex))
2215     return fMatchedTrackIndex->At(FindMatchedPosForCluster(clsIndex));
2216   else 
2217     return -1; 
2218 }
2219
2220 //__________________________________________________________
2221 Int_t AliEMCALRecoUtils::GetMatchedClusterIndex(Int_t trkIndex)
2222 {
2223   //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2224   //Get the index of matched cluster to this track
2225   //Works with ESDs and AODs
2226   
2227   if (IsTrackMatched(trkIndex))
2228     return fMatchedClusterIndex->At(FindMatchedPosForTrack(trkIndex));
2229   else 
2230     return -1; 
2231 }
2232
2233 //______________________________________________________________
2234 Bool_t AliEMCALRecoUtils::IsClusterMatched(Int_t clsIndex) const
2235 {
2236   //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2237   //Returns if the cluster has a match
2238   if (FindMatchedPosForCluster(clsIndex) < 999) 
2239     return kTRUE;
2240   else
2241     return kFALSE;
2242 }
2243
2244 //____________________________________________________________
2245 Bool_t AliEMCALRecoUtils::IsTrackMatched(Int_t trkIndex) const 
2246 {
2247   //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2248   //Returns if the track has a match
2249   if (FindMatchedPosForTrack(trkIndex) < 999) 
2250     return kTRUE;
2251   else
2252     return kFALSE;
2253 }
2254
2255 //______________________________________________________________________
2256 UInt_t AliEMCALRecoUtils::FindMatchedPosForCluster(Int_t clsIndex) const
2257 {
2258   //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2259   //Returns the position of the match in the fMatchedClusterIndex array
2260   Float_t tmpR = fCutR;
2261   UInt_t pos = 999;
2262   
2263   for (Int_t i=0; i<fMatchedClusterIndex->GetSize(); i++) 
2264   {
2265     if (fMatchedClusterIndex->At(i)==clsIndex) {
2266       Float_t r = TMath::Sqrt(fResidualEta->At(i)*fResidualEta->At(i) + fResidualPhi->At(i)*fResidualPhi->At(i));
2267       if (r<tmpR) {
2268         pos=i;
2269         tmpR=r;
2270         AliDebug(3,Form("Matched cluster index: index: %d, dEta: %2.4f, dPhi: %2.4f.\n",
2271                         fMatchedClusterIndex->At(i),fResidualEta->At(i),fResidualPhi->At(i)));
2272       }
2273     }
2274   }
2275   return pos;
2276 }
2277
2278 //____________________________________________________________________
2279 UInt_t AliEMCALRecoUtils::FindMatchedPosForTrack(Int_t trkIndex) const
2280 {
2281   //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2282   //Returns the position of the match in the fMatchedTrackIndex array
2283   Float_t tmpR = fCutR;
2284   UInt_t pos = 999;
2285   
2286   for (Int_t i=0; i<fMatchedTrackIndex->GetSize(); i++) 
2287   {
2288     if (fMatchedTrackIndex->At(i)==trkIndex) {
2289       Float_t r = TMath::Sqrt(fResidualEta->At(i)*fResidualEta->At(i) + fResidualPhi->At(i)*fResidualPhi->At(i));
2290       if (r<tmpR) {
2291         pos=i;
2292         tmpR=r;
2293         AliDebug(3,Form("Matched track index: index: %d, dEta: %2.4f, dPhi: %2.4f.\n",
2294                         fMatchedTrackIndex->At(i),fResidualEta->At(i),fResidualPhi->At(i)));
2295       }
2296     }
2297   }
2298   return pos;
2299 }
2300
2301 //__________________________________________________________________________
2302 Bool_t AliEMCALRecoUtils::IsGoodCluster(AliVCluster *cluster, 
2303                                         const AliEMCALGeometry *geom, 
2304                                         AliVCaloCells* cells, Int_t bc)
2305 {
2306   // check if the cluster survives some quality cut
2307   //
2308   //
2309   Bool_t isGood=kTRUE;
2310
2311   if (!cluster || !cluster->IsEMCAL())              return kFALSE;
2312   if (ClusterContainsBadChannel(geom,cluster->GetCellsAbsId(),cluster->GetNCells())) return kFALSE;
2313   if (!CheckCellFiducialRegion(geom,cluster,cells)) return kFALSE;
2314   if (IsExoticCluster(cluster, cells,bc))           return kFALSE;
2315
2316   return isGood;
2317 }
2318
2319 //__________________________________________________________
2320 Bool_t AliEMCALRecoUtils::IsAccepted(AliESDtrack *esdTrack)
2321 {
2322   // Given a esd track, return whether the track survive all the cuts
2323
2324   // The different quality parameter are first
2325   // retrieved from the track. then it is found out what cuts the
2326   // track did not survive and finally the cuts are imposed.
2327
2328   UInt_t status = esdTrack->GetStatus();
2329
2330   Int_t nClustersITS = esdTrack->GetITSclusters(0);
2331   Int_t nClustersTPC = esdTrack->GetTPCclusters(0);
2332
2333   Float_t chi2PerClusterITS = -1;
2334   Float_t chi2PerClusterTPC = -1;
2335   if (nClustersITS!=0)
2336     chi2PerClusterITS = esdTrack->GetITSchi2()/Float_t(nClustersITS);
2337   if (nClustersTPC!=0) 
2338     chi2PerClusterTPC = esdTrack->GetTPCchi2()/Float_t(nClustersTPC);
2339
2340
2341   //DCA cuts
2342   if (fTrackCutsType==kGlobalCut) {
2343     Float_t maxDCAToVertexXYPtDep = 0.0182 + 0.0350/TMath::Power(esdTrack->Pt(),1.01); //This expression comes from AliESDtrackCuts::GetStandardITSTPCTrackCuts2010()
2344     //AliDebug(3,Form("Track pT = %f, DCAtoVertexXY = %f",esdTrack->Pt(),MaxDCAToVertexXYPtDep));
2345     SetMaxDCAToVertexXY(maxDCAToVertexXYPtDep); //Set pT dependent DCA cut to vertex in x-y plane
2346   }
2347
2348   Float_t b[2];
2349   Float_t bCov[3];
2350   esdTrack->GetImpactParameters(b,bCov);
2351   if (bCov[0]<=0 || bCov[2]<=0) {
2352     AliDebug(1, "Estimated b resolution lower or equal zero!");
2353     bCov[0]=0; bCov[2]=0;
2354   }
2355
2356   Float_t dcaToVertexXY = b[0];
2357   Float_t dcaToVertexZ = b[1];
2358   Float_t dcaToVertex = -1;
2359
2360   if (fCutDCAToVertex2D)
2361     dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY/fCutMaxDCAToVertexXY/fCutMaxDCAToVertexXY + dcaToVertexZ*dcaToVertexZ/fCutMaxDCAToVertexZ/fCutMaxDCAToVertexZ);
2362   else
2363     dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY + dcaToVertexZ*dcaToVertexZ);
2364     
2365   // cut the track?
2366   
2367   Bool_t cuts[kNCuts];
2368   for (Int_t i=0; i<kNCuts; i++) cuts[i]=kFALSE;
2369   
2370   // track quality cuts
2371   if (fCutRequireTPCRefit && (status&AliESDtrack::kTPCrefit)==0)
2372     cuts[0]=kTRUE;
2373   if (fCutRequireITSRefit && (status&AliESDtrack::kITSrefit)==0)
2374     cuts[1]=kTRUE;
2375   if (nClustersTPC<fCutMinNClusterTPC)
2376     cuts[2]=kTRUE;
2377   if (nClustersITS<fCutMinNClusterITS) 
2378     cuts[3]=kTRUE;
2379   if (chi2PerClusterTPC>fCutMaxChi2PerClusterTPC) 
2380     cuts[4]=kTRUE; 
2381   if (chi2PerClusterITS>fCutMaxChi2PerClusterITS) 
2382     cuts[5]=kTRUE;  
2383   if (!fCutAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0)
2384     cuts[6]=kTRUE;
2385   if (fCutDCAToVertex2D && dcaToVertex > 1)
2386     cuts[7] = kTRUE;
2387   if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexXY) > fCutMaxDCAToVertexXY)
2388     cuts[8] = kTRUE;
2389   if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexZ) > fCutMaxDCAToVertexZ)
2390     cuts[9] = kTRUE;
2391
2392   if (fTrackCutsType==kGlobalCut) {
2393     //Require at least one SPD point + anything else in ITS
2394     if ( (esdTrack->HasPointOnITSLayer(0) || esdTrack->HasPointOnITSLayer(1)) == kFALSE)
2395       cuts[10] = kTRUE;
2396   }
2397
2398   // ITS
2399   if (fCutRequireITSStandAlone || fCutRequireITSpureSA) {
2400     if ((status & AliESDtrack::kITSin) == 0 || (status & AliESDtrack::kTPCin)) {
2401       // TPC tracks
2402       cuts[11] = kTRUE; 
2403     } else {
2404       // ITS standalone tracks
2405       if (fCutRequireITSStandAlone && !fCutRequireITSpureSA) {
2406         if (status & AliESDtrack::kITSpureSA) cuts[11] = kTRUE;
2407       } else if (fCutRequireITSpureSA) {
2408         if (!(status & AliESDtrack::kITSpureSA)) cuts[11] = kTRUE;
2409       }
2410     }
2411   }
2412   
2413   Bool_t cut=kFALSE;
2414   for (Int_t i=0; i<kNCuts; i++)
2415     if (cuts[i]) { cut = kTRUE ; }
2416
2417     // cut the track
2418   if (cut) 
2419     return kFALSE;
2420   else 
2421     return kTRUE;
2422 }
2423
2424 //_____________________________________
2425 void AliEMCALRecoUtils::InitTrackCuts()
2426 {
2427   //Intilize the track cut criteria
2428   //By default these cuts are set according to AliESDtrackCuts::GetStandardTPCOnlyTrackCuts()
2429   //Also you can customize the cuts using the setters
2430   
2431   switch (fTrackCutsType)
2432   {
2433     case kTPCOnlyCut:
2434     {
2435       AliInfo(Form("Track cuts for matching: GetStandardTPCOnlyTrackCuts()"));
2436       //TPC
2437       SetMinNClustersTPC(70);
2438       SetMaxChi2PerClusterTPC(4);
2439       SetAcceptKinkDaughters(kFALSE);
2440       SetRequireTPCRefit(kFALSE);
2441       
2442       //ITS
2443       SetRequireITSRefit(kFALSE);
2444       SetMaxDCAToVertexZ(3.2);
2445       SetMaxDCAToVertexXY(2.4);
2446       SetDCAToVertex2D(kTRUE);
2447       
2448       break;
2449     }
2450       
2451     case kGlobalCut:
2452     {
2453       AliInfo(Form("Track cuts for matching: GetStandardITSTPCTrackCuts2010(kTURE)"));
2454       //TPC
2455       SetMinNClustersTPC(70);
2456       SetMaxChi2PerClusterTPC(4);
2457       SetAcceptKinkDaughters(kFALSE);
2458       SetRequireTPCRefit(kTRUE);
2459       
2460       //ITS
2461       SetRequireITSRefit(kTRUE);
2462       SetMaxDCAToVertexZ(2);
2463       SetMaxDCAToVertexXY();
2464       SetDCAToVertex2D(kFALSE);
2465       
2466       break;
2467     }
2468       
2469     case kLooseCut:
2470     {
2471       AliInfo(Form("Track cuts for matching: Loose cut w/o DCA cut"));
2472       SetMinNClustersTPC(50);
2473       SetAcceptKinkDaughters(kTRUE);
2474       
2475       break;
2476     }
2477
2478     case kITSStandAlone:
2479     {
2480       AliInfo(Form("Track cuts for matching: ITS Stand Alone tracks cut w/o DCA cut"));
2481       SetRequireITSRefit(kTRUE);
2482       SetRequireITSStandAlone(kTRUE);
2483       SetITSTrackSA(kTRUE);
2484       break;
2485     }
2486     
2487   }
2488 }
2489
2490
2491 //________________________________________________________________________
2492 void AliEMCALRecoUtils::SetClusterMatchedToTrack(const AliVEvent *event)
2493 {
2494   // Checks if tracks are matched to EMC clusters and set the matched EMCAL cluster index to ESD track. 
2495
2496   Int_t nTracks = event->GetNumberOfTracks();
2497   for (Int_t iTrack = 0; iTrack < nTracks; ++iTrack) 
2498   {
2499     AliVTrack* track = dynamic_cast<AliVTrack*>(event->GetTrack(iTrack));
2500     if (!track) 
2501     {
2502       AliWarning(Form("Could not receive track %d", iTrack));
2503       continue;
2504     }
2505     
2506     Int_t matchClusIndex = GetMatchedClusterIndex(iTrack);       
2507     track->SetEMCALcluster(matchClusIndex); //sets -1 if track not matched within residual
2508     /*the following can be done better if AliVTrack::SetStatus will be there. Patch pending with Andreas/Peter*/
2509     AliESDtrack* esdtrack = dynamic_cast<AliESDtrack*>(track);
2510     if (esdtrack) { 
2511       if (matchClusIndex != -1) 
2512         esdtrack->SetStatus(AliESDtrack::kEMCALmatch);
2513       else
2514         esdtrack->ResetStatus(AliESDtrack::kEMCALmatch);
2515     } else {
2516       AliAODTrack* aodtrack = dynamic_cast<AliAODTrack*>(track);
2517       if (matchClusIndex != -1) 
2518         aodtrack->SetStatus(AliESDtrack::kEMCALmatch);
2519       else
2520         aodtrack->ResetStatus(AliESDtrack::kEMCALmatch);
2521     }
2522   }
2523   AliDebug(2,"Track matched to closest cluster");  
2524 }
2525
2526 //_________________________________________________________________________
2527 void AliEMCALRecoUtils::SetTracksMatchedToCluster(const AliVEvent *event)
2528 {
2529   // Checks if EMC clusters are matched to ESD track.
2530   // Adds track indexes of all the tracks matched to a cluster withing residuals in ESDCalocluster.
2531   
2532   for (Int_t iClus=0; iClus < event->GetNumberOfCaloClusters(); ++iClus) 
2533   {
2534     AliVCluster *cluster = event->GetCaloCluster(iClus);
2535     if (!cluster->IsEMCAL()) 
2536       continue;
2537     
2538     Int_t nTracks = event->GetNumberOfTracks();
2539     TArrayI arrayTrackMatched(nTracks);
2540     
2541     // Get the closest track matched to the cluster
2542     Int_t nMatched = 0;
2543     Int_t matchTrackIndex = GetMatchedTrackIndex(iClus);
2544     if (matchTrackIndex != -1) 
2545     {
2546       arrayTrackMatched[nMatched] = matchTrackIndex;
2547       nMatched++;
2548     }
2549     
2550     // Get all other tracks matched to the cluster
2551     for (Int_t iTrk=0; iTrk<nTracks; ++iTrk) 
2552     {
2553       AliVTrack* track = dynamic_cast<AliVTrack*>(event->GetTrack(iTrk));
2554       
2555       if( !track ) continue;
2556       
2557       if ( iTrk == matchTrackIndex ) continue;
2558       
2559       if ( track->GetEMCALcluster() == iClus )
2560       {
2561         arrayTrackMatched[nMatched] = iTrk;
2562         ++nMatched;
2563       }
2564     }
2565     
2566     //printf("Tender::SetTracksMatchedToCluster - cluster E %f, N matches %d, first match %d\n",cluster->E(),nMatched,arrayTrackMatched[0]);
2567     
2568     arrayTrackMatched.Set(nMatched);
2569     AliESDCaloCluster *esdcluster = dynamic_cast<AliESDCaloCluster*>(cluster);
2570     if (esdcluster) 
2571       esdcluster->AddTracksMatched(arrayTrackMatched);
2572     else if (nMatched>0) {
2573       AliAODCaloCluster *aodcluster = dynamic_cast<AliAODCaloCluster*>(cluster);
2574       if (aodcluster)
2575         aodcluster->AddTrackMatched(event->GetTrack(arrayTrackMatched.At(0)));
2576     }
2577     
2578     Float_t eta= -999, phi = -999;
2579     if (matchTrackIndex != -1) 
2580       GetMatchedResiduals(iClus, eta, phi);
2581     cluster->SetTrackDistance(phi, eta);
2582   }
2583   
2584   AliDebug(2,"Cluster matched to tracks");  
2585 }
2586
2587 //___________________________________________________
2588 void AliEMCALRecoUtils::Print(const Option_t *) const 
2589 {
2590   // Print Parameters
2591   
2592   printf("AliEMCALRecoUtils Settings: \n");
2593   printf("Misalignment shifts\n");
2594   for (Int_t i=0; i<5; i++) printf("\t sector %d, traslation (x,y,z)=(%f,%f,%f), rotation (x,y,z)=(%f,%f,%f)\n",i, 
2595                                   fMisalTransShift[i*3],fMisalTransShift[i*3+1],fMisalTransShift[i*3+2],
2596                                   fMisalRotShift[i*3],  fMisalRotShift[i*3+1],  fMisalRotShift[i*3+2]   );
2597   printf("Non linearity function %d, parameters:\n", fNonLinearityFunction);
2598   for (Int_t i=0; i<6; i++) printf("param[%d]=%f\n",i, fNonLinearityParams[i]);
2599   
2600   printf("Position Recalculation option %d, Particle Type %d, fW0 %2.2f, Recalibrate Data %d \n",fPosAlgo,fParticleType,fW0, fRecalibration);
2601
2602   printf("Matching criteria: ");
2603   if (fCutEtaPhiSum) {
2604     printf("sqrt(dEta^2+dPhi^2)<%4.3f\n",fCutR);
2605   } else if (fCutEtaPhiSeparate) {
2606     printf("dEta<%4.3f, dPhi<%4.3f\n",fCutEta,fCutPhi);
2607   } else {
2608     printf("Error\n");
2609     printf("please specify your cut criteria\n");
2610     printf("To cut on sqrt(dEta^2+dPhi^2), use: SwitchOnCutEtaPhiSum()\n");
2611     printf("To cut on dEta and dPhi separately, use: SwitchOnCutEtaPhiSeparate()\n");
2612   }
2613
2614   printf("Mass hypothesis = %2.3f [GeV/c^2], extrapolation step to surface = %2.2f[cm], step to cluster = %2.2f[cm]\n",fMass,fStepSurface, fStepCluster);
2615   printf("Cluster selection window: dR < %2.0f\n",fClusterWindow);
2616
2617   printf("Track cuts: \n");
2618   printf("Minimum track pT: %1.2f\n",fCutMinTrackPt);
2619   printf("AOD track selection: tpc only %d, or hybrid %d, or mask: %d\n",fAODTPCOnlyTracks,fAODHybridTracks, fAODFilterMask);
2620   printf("TPCRefit = %d, ITSRefit = %d\n",fCutRequireTPCRefit,fCutRequireITSRefit);
2621   printf("AcceptKinks = %d\n",fCutAcceptKinkDaughters);
2622   printf("MinNCulsterTPC = %d, MinNClusterITS = %d\n",fCutMinNClusterTPC,fCutMinNClusterITS);
2623   printf("MaxChi2TPC = %2.2f, MaxChi2ITS = %2.2f\n",fCutMaxChi2PerClusterTPC,fCutMaxChi2PerClusterITS);
2624   printf("DCSToVertex2D = %d, MaxDCAToVertexXY = %2.2f, MaxDCAToVertexZ = %2.2f\n",fCutDCAToVertex2D,fCutMaxDCAToVertexXY,fCutMaxDCAToVertexZ);
2625 }