Move exotic cell rejection methods from AliAnalysisTaskEMCALClusterize to AliEMCALRec...
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALRecoUtils.h
1 #ifndef ALIEMCALRECOUTILS_H
2 #define ALIEMCALRECOUTILS_H
3
4 /* $Id: AliEMCALRecoUtils.h 33808 2009-07-15 09:48:08Z gconesab $ */
5
6 ///////////////////////////////////////////////////////////////////////////////
7 //
8 // Class AliEMCALRecoUtils
9 // Some utilities to recalculate the cluster position or energy linearity
10 //
11 //
12 // Author:  Gustavo Conesa (LPSC- Grenoble) 
13 //          Track matching part: Rongrong Ma (Yale)
14 ///////////////////////////////////////////////////////////////////////////////
15
16 //Root includes
17 #include <TNamed.h>
18 #include <TMath.h>
19 class TObjArray;
20 class TArrayI;
21 class TArrayF;
22 #include <TH2I.h>
23 class TH2F;
24 #include <TRandom3.h>
25
26 //AliRoot includes
27 class AliVCluster;
28 class AliVCaloCells;
29 class AliVEvent;
30
31 // EMCAL includes
32 class AliEMCALGeometry;
33 class AliEMCALPIDUtils;
34 class AliESDtrack;
35 class AliExternalTrackParam;
36
37 class AliEMCALRecoUtils : public TNamed {
38   
39 public:
40   
41   AliEMCALRecoUtils();
42   AliEMCALRecoUtils(const AliEMCALRecoUtils&); 
43   AliEMCALRecoUtils& operator=(const AliEMCALRecoUtils&); 
44   virtual ~AliEMCALRecoUtils() ;  
45   void     Print(const Option_t*) const;
46
47   //enums
48   enum     NonlinearityFunctions{kPi0MC=0,kPi0GammaGamma=1,kPi0GammaConversion=2,kNoCorrection=3,kBeamTest=4,kBeamTestCorrected=5};
49   enum     PositionAlgorithms{kUnchanged=-1,kPosTowerIndex=0, kPosTowerGlobal=1};
50   enum     ParticleType{kPhoton=0, kElectron=1,kHadron =2, kUnknown=-1};
51   enum     { kNCuts = 11 }; //track matching
52   enum     TrackCutsType{kTPCOnlyCut=0, kGlobalCut=1, kLooseCut=2};
53
54   //-----------------------------------------------------
55   //Position recalculation
56   //-----------------------------------------------------
57
58   void     RecalculateClusterPosition               (AliEMCALGeometry *geom, AliVCaloCells* cells, AliVCluster* clu); 
59   void     RecalculateClusterPositionFromTowerIndex (AliEMCALGeometry *geom, AliVCaloCells* cells, AliVCluster* clu); 
60   void     RecalculateClusterPositionFromTowerGlobal(AliEMCALGeometry *geom, AliVCaloCells* cells, AliVCluster* clu); 
61   
62   Float_t  GetCellWeight(const Float_t eCell, const Float_t eCluster) const { return TMath::Max( 0., fW0 + TMath::Log( eCell / eCluster )) ; }
63   
64   Float_t  GetDepth(const Float_t eCluster, const Int_t iParticle, const Int_t iSM) const ; 
65   
66   void     GetMaxEnergyCell(AliEMCALGeometry *geom, AliVCaloCells* cells, AliVCluster* clu, 
67                             Int_t & absId,  Int_t& iSupMod, Int_t& ieta, Int_t& iphi, Bool_t &shared);
68   
69   Float_t  GetMisalTransShift(const Int_t i)       const { if(i < 15 ) { return fMisalTransShift[i] ; }
70                                                            else        { AliInfo(Form("Index %d larger than 15, do nothing\n",i)) ; 
71                                                                          return 0.                  ; } }
72   Float_t* GetMisalTransShiftArray()                     { return fMisalTransShift ; }
73
74   void     SetMisalTransShift(const Int_t i, const Float_t shift) {
75                                                            if(i < 15 ) { fMisalTransShift[i] = shift ; }
76                                                            else        { AliInfo(Form("Index %d larger than 15, do nothing\n",i)) ; } }
77   void     SetMisalTransShiftArray(Float_t * misal)               { for(Int_t i = 0; i < 15; i++) fMisalTransShift[i] = misal[i]  ; }
78
79   Float_t  GetMisalRotShift(const Int_t i)         const { if(i < 15 ) { return fMisalRotShift[i]    ; }
80                                                            else        { AliInfo(Form("Index %d larger than 15, do nothing\n",i)) ; 
81                                                                          return 0.                   ; } }
82   
83   Float_t* GetMisalRotShiftArray()                       { return fMisalRotShift                     ; }
84   
85   void     SetMisalRotShift(const Int_t i, const Float_t shift) {
86                                                            if(i < 15 ) { fMisalRotShift[i] = shift   ; }
87                                                            else        { AliInfo(Form("Index %d larger than 15, do nothing\n",i)) ; } }
88   
89   void     SetMisalRotShiftArray(Float_t * misal)        { for(Int_t i = 0; i < 15; i++)fMisalRotShift[i] = misal[i] ; }
90   
91   Int_t    GetParticleType()                       const { return  fParticleType    ; }
92   void     SetParticleType(Int_t particle)               { fParticleType = particle ; }
93   
94   Int_t    GetPositionAlgorithm()                  const { return fPosAlgo          ; }
95   void     SetPositionAlgorithm(Int_t alg)               { fPosAlgo = alg           ; }
96   
97   Float_t  GetW0()                                 const { return fW0               ; }
98   void     SetW0(Float_t w0)                             { fW0  = w0                ; }
99
100   //-----------------------------------------------------
101   // Non Linearity
102   //-----------------------------------------------------
103
104   Float_t  CorrectClusterEnergyLinearity(AliVCluster* clu) ;
105   
106   Float_t  GetNonLinearityParam(const Int_t i)     const { if(i < 7 ){ return fNonLinearityParams[i] ; }
107                                                           else      { AliInfo(Form("Index %d larger than 7, do nothing\n",i)) ; 
108                                                                        return 0.                     ; } }
109   void     SetNonLinearityParam(const Int_t i, const Float_t param) {
110                                                           if(i < 7 ){fNonLinearityParams[i] = param ; }
111                                                           else { AliInfo(Form("Index %d larger than 7, do nothing\n",i)) ; } }
112   void     InitNonLinearityParam();
113
114   Int_t    GetNonLinearityFunction() const               { return fNonLinearityFunction    ; }
115   void     SetNonLinearityFunction(Int_t fun)            { fNonLinearityFunction = fun     ; InitNonLinearityParam() ; }
116
117   void     SetNonLinearityThreshold(Int_t threshold)     { fNonLinearThreshold = threshold ; } //only for Alexie's non linearity correction
118   Int_t    GetNonLinearityThreshold()              const { return fNonLinearThreshold      ; }
119 //  
120   //-----------------------------------------------------
121   // MC clusters energy smearing
122   //-----------------------------------------------------
123   
124   Float_t  SmearClusterEnergy(AliVCluster* clu) ;
125   void     SwitchOnClusterEnergySmearing()               { fSmearClusterEnergy = kTRUE         ; }
126   void     SwitchOffClusterEnergySmearing()              { fSmearClusterEnergy = kFALSE        ; }
127   Bool_t   IsClusterEnergySmeared()                const { return fSmearClusterEnergy          ; }   
128   void     SetSmearingParameters(Int_t i, Float_t param) { if(i < 3){ fSmearClusterParam[i] = param ; }
129                                                            else     { AliInfo(Form("Index %d larger than 2, do nothing\n",i)) ; } }
130   //-----------------------------------------------------
131   // Recalibration
132   //-----------------------------------------------------
133   Bool_t   AcceptCalibrateCell(const Int_t absId, const Int_t bc,
134                                Float_t & amp, Double_t & time, AliVCaloCells* cells) ; // Energy and Time
135   void     RecalibrateCells(AliVCaloCells * cells, Int_t bc) ; // Energy and Time
136   void     RecalibrateClusterEnergy(AliEMCALGeometry* geom, AliVCluster* cluster, AliVCaloCells * cells, const Int_t bc=-1) ; // Energy and time
137
138   // Energy recalibration
139   Bool_t   IsRecalibrationOn()                     const { return fRecalibration ; }
140   void     SwitchOffRecalibration()                      { fRecalibration = kFALSE ; }
141   void     SwitchOnRecalibration()                       { fRecalibration = kTRUE  ; 
142                                                            if(!fEMCALRecalibrationFactors)InitEMCALRecalibrationFactors() ; }
143   void     InitEMCALRecalibrationFactors() ;
144
145   TH2F *   GetEMCALChannelRecalibrationFactors(Int_t iSM)     const { return (TH2F*)fEMCALRecalibrationFactors->At(iSM) ; }     
146   void     SetEMCALChannelRecalibrationFactors(TObjArray *map)      { fEMCALRecalibrationFactors = map                  ; }
147   void     SetEMCALChannelRecalibrationFactors(Int_t iSM , TH2F* h) { fEMCALRecalibrationFactors->AddAt(h,iSM)          ; }
148   
149   Float_t  GetEMCALChannelRecalibrationFactor(Int_t iSM , Int_t iCol, Int_t iRow) const { 
150     if(fEMCALRecalibrationFactors) 
151       return (Float_t) ((TH2F*)fEMCALRecalibrationFactors->At(iSM))->GetBinContent(iCol,iRow); 
152     else return 1 ; } 
153         
154   void     SetEMCALChannelRecalibrationFactor(Int_t iSM , Int_t iCol, Int_t iRow, Double_t c = 1) { 
155     if(!fEMCALRecalibrationFactors) InitEMCALRecalibrationFactors() ;
156     ((TH2F*)fEMCALRecalibrationFactors->At(iSM))->SetBinContent(iCol,iRow,c) ; }
157   
158   //Recalibrate channels energy with run dependent corrections
159   void     SwitchOffRunDepCorrection()                   { fUseRunCorrectionFactors = kFALSE ; }
160   void     SwitchOnRunDepCorrection()                    { fUseRunCorrectionFactors = kTRUE  ; 
161                                                            SwitchOnRecalibration()           ; }
162   void     SetRunDependentCorrections(Int_t runnumber);
163       
164   // Time Recalibration  
165   void     RecalibrateCellTime(const Int_t absId, const Int_t bc, Double_t & time);
166   
167   Bool_t   IsTimeRecalibrationOn()                 const { return fTimeRecalibration   ; }
168   void     SwitchOffTimeRecalibration()                  { fTimeRecalibration = kFALSE ; }
169   void     SwitchOnTimeRecalibration()                   { fTimeRecalibration = kTRUE  ; 
170     if(!fEMCALTimeRecalibrationFactors)InitEMCALTimeRecalibrationFactors() ; }
171   void     InitEMCALTimeRecalibrationFactors() ;
172   
173   Float_t  GetEMCALChannelTimeRecalibrationFactor(const Int_t bc, const Int_t absID) const { 
174     if(fEMCALTimeRecalibrationFactors) 
175       return (Float_t) ((TH1F*)fEMCALTimeRecalibrationFactors->At(bc))->GetBinContent(absID); 
176     else return 0 ; } 
177         
178   void     SetEMCALChannelTimeRecalibrationFactor(const Int_t bc, const Int_t absID, Double_t c = 0) { 
179     if(!fEMCALTimeRecalibrationFactors) InitEMCALTimeRecalibrationFactors() ;
180     ((TH1F*)fEMCALTimeRecalibrationFactors->At(bc))->SetBinContent(absID,c) ; }  
181   
182   TH1F *   GetEMCALChannelTimeRecalibrationFactors(const Int_t bc)const       { return (TH1F*)fEMCALTimeRecalibrationFactors->At(bc) ; }        
183   void     SetEMCALChannelTimeRecalibrationFactors(TObjArray *map)            { fEMCALTimeRecalibrationFactors = map                 ; }
184   void     SetEMCALChannelTimeRecalibrationFactors(const Int_t bc , TH1F* h)  { fEMCALTimeRecalibrationFactors->AddAt(h,bc)          ; }
185   
186   //-----------------------------------------------------
187   // Modules fiducial region, remove clusters in borders
188   //-----------------------------------------------------
189
190   Bool_t   CheckCellFiducialRegion(AliEMCALGeometry* geom, AliVCluster* cluster, AliVCaloCells* cells) ;
191   void     SetNumberOfCellsFromEMCALBorder(const Int_t n){ fNCellsFromEMCALBorder = n      ; }
192   Int_t    GetNumberOfCellsFromEMCALBorder()      const  { return fNCellsFromEMCALBorder   ; }
193     
194   void     SwitchOnNoFiducialBorderInEMCALEta0()         { fNoEMCALBorderAtEta0 = kTRUE    ; }
195   void     SwitchOffNoFiducialBorderInEMCALEta0()        { fNoEMCALBorderAtEta0 = kFALSE   ; }
196   Bool_t   IsEMCALNoBorderAtEta0()                 const { return fNoEMCALBorderAtEta0     ; }
197   
198   //-----------------------------------------------------
199   // Bad channels
200   //-----------------------------------------------------
201
202   Bool_t   IsBadChannelsRemovalSwitchedOn()        const { return fRemoveBadChannels       ; }
203   void     SwitchOffBadChannelsRemoval()                 { fRemoveBadChannels = kFALSE     ; }
204   void     SwitchOnBadChannelsRemoval ()                 { fRemoveBadChannels = kTRUE ; 
205                                                            if(!fEMCALBadChannelMap)InitEMCALBadChannelStatusMap() ; }
206         
207   Bool_t   IsDistanceToBadChannelRecalculated()    const { return fRecalDistToBadChannels   ; }
208   void     SwitchOffDistToBadChannelRecalculation()      { fRecalDistToBadChannels = kFALSE ; }
209   void     SwitchOnDistToBadChannelRecalculation()       { fRecalDistToBadChannels = kTRUE  ; 
210                                                            if(!fEMCALBadChannelMap)InitEMCALBadChannelStatusMap() ; }
211   
212   void     InitEMCALBadChannelStatusMap() ;
213         
214   Int_t    GetEMCALChannelStatus(Int_t iSM , Int_t iCol, Int_t iRow) const { 
215     if(fEMCALBadChannelMap) return (Int_t) ((TH2I*)fEMCALBadChannelMap->At(iSM))->GetBinContent(iCol,iRow); 
216     else return 0;}//Channel is ok by default
217         
218   void     SetEMCALChannelStatus(Int_t iSM , Int_t iCol, Int_t iRow, Double_t c = 1) { 
219                                                            if(!fEMCALBadChannelMap)InitEMCALBadChannelStatusMap()               ;
220                                                            ((TH2I*)fEMCALBadChannelMap->At(iSM))->SetBinContent(iCol,iRow,c)    ; }
221         
222   TH2I *   GetEMCALChannelStatusMap(Int_t iSM)     const { return (TH2I*)fEMCALBadChannelMap->At(iSM) ; }
223   void     SetEMCALChannelStatusMap(TObjArray *map)      { fEMCALBadChannelMap = map                  ; }
224   void     SetEMCALChannelStatusMap(Int_t iSM , TH2I* h) { fEMCALBadChannelMap->AddAt(h,iSM)          ; }
225
226   Bool_t   ClusterContainsBadChannel(AliEMCALGeometry* geom, UShort_t* cellList, const Int_t nCells);
227  
228   //-----------------------------------------------------
229   // Recalculate other cluster parameters
230   //-----------------------------------------------------
231
232   void     RecalculateClusterDistanceToBadChannel (AliEMCALGeometry * geom, AliVCaloCells* cells, AliVCluster * cluster);
233   void     RecalculateClusterShowerShapeParameters(AliEMCALGeometry * geom, AliVCaloCells* cells, AliVCluster * cluster);
234   void     RecalculateClusterPID(AliVCluster * cluster);
235
236   AliEMCALPIDUtils * GetPIDUtils() { return fPIDUtils;}
237
238
239   //----------------------------------------------------
240   // Track matching
241   //----------------------------------------------------
242
243   void     FindMatches(AliVEvent *event, TObjArray * clusterArr=0x0, AliEMCALGeometry *geom=0x0);
244   Int_t    FindMatchedClusterInEvent(AliESDtrack *track, AliVEvent *event, AliEMCALGeometry *geom, Float_t &dEta, Float_t &dPhi);
245   Int_t    FindMatchedClusterInClusterArr(AliExternalTrackParam *emcalParam, AliExternalTrackParam *trkParam, TObjArray * clusterArr, Float_t &dEta, Float_t &dPhi);
246   
247   static Bool_t ExtrapolateTrackToEMCalSurface(AliExternalTrackParam *trkParam, Double_t emcalR, Double_t mass, Double_t step, Float_t &eta, Float_t &phi);
248   static Bool_t ExtrapolateTrackToPosition(AliExternalTrackParam *trkParam, Float_t *clsPos, Double_t mass, Double_t step, Float_t &tmpEta, Float_t &tmpPhi);
249   static Bool_t ExtrapolateTrackToCluster(AliExternalTrackParam *trkParam, AliVCluster *cluster, Double_t mass, Double_t step, Float_t &tmpEta, Float_t &tmpPhi);
250   Bool_t ExtrapolateTrackToCluster(AliExternalTrackParam *trkParam, AliVCluster *cluster, Float_t &tmpEta, Float_t &tmpPhi);
251
252   UInt_t   FindMatchedPosForCluster(Int_t clsIndex) const;
253   UInt_t   FindMatchedPosForTrack(Int_t trkIndex)   const;
254   
255   void     GetMatchedResiduals(Int_t clsIndex, Float_t &dEta, Float_t &dPhi);
256   void     GetMatchedClusterResiduals(Int_t trkIndex, Float_t &dEta, Float_t &dPhi);
257   Int_t    GetMatchedTrackIndex(Int_t clsIndex);
258   Int_t    GetMatchedClusterIndex(Int_t trkIndex);
259   
260   Bool_t   IsClusterMatched(Int_t clsIndex)         const;
261   Bool_t   IsTrackMatched(Int_t trkIndex)           const;
262
263
264   void     SwitchOnCutEtaPhiSum()                     { fCutEtaPhiSum      = kTRUE    ; 
265                                                         fCutEtaPhiSeparate = kFALSE   ; }
266   void     SwitchOnCutEtaPhiSeparate()                { fCutEtaPhiSeparate = kTRUE    ;
267                                                         fCutEtaPhiSum      = kFALSE   ; }
268
269   Float_t  GetCutR()                            const { return fCutR                  ; }
270   Float_t  GetCutEta()                          const { return fCutEta                ; }
271   Float_t  GetCutPhi()                          const { return fCutPhi                ; }
272   Double_t GetClusterWindow()                   const { return fClusterWindow         ; }
273   void     SetCutR(Float_t cutR)                      { fCutR   = cutR                ; }
274   void     SetCutEta(Float_t cutEta)                  { fCutEta = cutEta              ; }
275   void     SetCutPhi(Float_t cutPhi)                  { fCutPhi = cutPhi              ; }
276   void     SetClusterWindow(Double_t window)          { fClusterWindow = window       ; }
277   void     SetCutZ(Float_t cutZ)                      { printf("Obsolete fucntion of cutZ=%1.1f\n",cutZ) ; } //Obsolete
278
279   Double_t GetMass()                            const { return fMass                  ; }
280   Double_t GetStep()                            const { return fStepCluster           ; }
281   Double_t GetStepSurface()                     const { return fStepSurface           ; }
282   void     SetMass(Double_t mass)                     { fMass = mass                  ; }
283   void     SetStep(Double_t step)                     { fStepCluster = step           ; }
284   void     SetStepSurface(Double_t step)              { fStepSurface = step           ; }
285  
286   // Exotic cells / clusters
287   
288   Bool_t   IsExoticCell(const Int_t absId, AliVCaloCells* cells, const Int_t bc =-1) ;
289   void     SwitchOnRejectExoticCell()                 { fRejectExoticCells = kTRUE     ; }
290   void     SwitchOffRejectExoticCell()                { fRejectExoticCells = kFALSE    ; } 
291    
292   void     SetExoticCellFractionCut(Float_t f)        { fExoticCellFraction     = f    ; }
293   void     SetExoticCellDiffTimeCut(Float_t dt)       { fExoticCellDiffTime     = dt   ; }
294   void     SetExoticCellMinAmplitudeCut(Float_t ma)   { fExoticCellMinAmplitude = ma   ; }
295   
296   Bool_t   IsExoticCluster(AliVCluster *cluster, AliVCaloCells* cells, const Int_t bc=0) ;
297   void     SwitchOnRejectExoticCluster()              { fRejectExoticCluster = kTRUE   ;
298                                                         fRejectExoticCells   = kTRUE   ; }
299   void     SwitchOffRejectExoticCluster()             { fRejectExoticCluster = kFALSE  ; }
300   Bool_t   IsRejectExoticCluster()              const { return fRejectExoticCluster    ; }
301   
302   //Cluster cut
303   Bool_t   IsGoodCluster(AliVCluster *cluster, AliEMCALGeometry *geom, AliVCaloCells* cells, const Int_t bc =-1);
304
305   //Track Cuts 
306   Bool_t   IsAccepted(AliESDtrack *track);
307   void     InitTrackCuts();
308   void     SetTrackCutsType(Int_t type)              { fTrackCutsType = type           ; 
309                                                        InitTrackCuts()                 ; }
310   Int_t    GetTrackCutsType() const                  { return fTrackCutsType; }
311
312   // track quality cut setters  
313   void     SetMinTrackPt(Double_t pt=0)              { fCutMinTrackPt           = pt   ; }
314   void     SetMinNClustersTPC(Int_t min=-1)          { fCutMinNClusterTPC       = min  ; }
315   void     SetMinNClustersITS(Int_t min=-1)          { fCutMinNClusterITS       = min  ; }
316   void     SetMaxChi2PerClusterTPC(Float_t max=1e10) { fCutMaxChi2PerClusterTPC = max  ; }
317   void     SetMaxChi2PerClusterITS(Float_t max=1e10) { fCutMaxChi2PerClusterITS = max  ; }
318   void     SetRequireTPCRefit(Bool_t b=kFALSE)       { fCutRequireTPCRefit      = b    ; }
319   void     SetRequireITSRefit(Bool_t b=kFALSE)       { fCutRequireITSRefit      = b    ; }
320   void     SetAcceptKinkDaughters(Bool_t b=kTRUE)    { fCutAcceptKinkDaughters  = b    ; }
321   void     SetMaxDCAToVertexXY(Float_t dist=1e10)    { fCutMaxDCAToVertexXY     = dist ; }
322   void     SetMaxDCAToVertexZ(Float_t dist=1e10)     { fCutMaxDCAToVertexZ      = dist ; }
323   void     SetDCAToVertex2D(Bool_t b=kFALSE)         { fCutDCAToVertex2D        = b    ; }
324
325   // getters                                                            
326   Double_t GetMinTrackPt()                     const { return fCutMinTrackPt           ; }
327   Int_t    GetMinNClusterTPC()                 const { return fCutMinNClusterTPC       ; }
328   Int_t    GetMinNClustersITS()                const { return fCutMinNClusterITS       ; }
329   Float_t  GetMaxChi2PerClusterTPC()           const { return fCutMaxChi2PerClusterTPC ; }
330   Float_t  GetMaxChi2PerClusterITS()           const { return fCutMaxChi2PerClusterITS ; }
331   Bool_t   GetRequireTPCRefit()                const { return fCutRequireTPCRefit      ; }
332   Bool_t   GetRequireITSRefit()                const { return fCutRequireITSRefit      ; }
333   Bool_t   GetAcceptKinkDaughters()            const { return fCutAcceptKinkDaughters  ; }
334   Float_t  GetMaxDCAToVertexXY()               const { return fCutMaxDCAToVertexXY     ; }
335   Float_t  GetMaxDCAToVertexZ()                const { return fCutMaxDCAToVertexZ      ; }
336   Bool_t   GetDCAToVertex2D()                  const { return fCutDCAToVertex2D        ; }
337
338
339 private:  
340   //Position recalculation
341   Float_t    fMisalTransShift[15];       // Shift parameters
342   Float_t    fMisalRotShift[15];         // Shift parameters
343   Int_t      fParticleType;              // Particle type for depth calculation
344   Int_t      fPosAlgo;                   // Position recalculation algorithm
345   Float_t    fW0;                        // Weight0
346     
347   // Non linearity
348   Int_t      fNonLinearityFunction;      // Non linearity function choice
349   Float_t    fNonLinearityParams[7];     // Parameters for the non linearity function
350   Int_t      fNonLinearThreshold;        // Non linearity threshold value for kBeamTesh non linearity function 
351   
352   // Energy smearing for MC
353   Bool_t     fSmearClusterEnergy;        // Smear cluster energy, to be done only for simulated data to match real data
354   Float_t    fSmearClusterParam[3];      // Smearing parameters
355   TRandom3   fRandom;                    // Random generator
356     
357   // Energy Recalibration 
358   Bool_t     fCellsRecalibrated;         // Internal bool to check if cells (time/energy) where recalibrated and not recalibrate them when recalculating different things
359   Bool_t     fRecalibration;             // Switch on or off the recalibration
360   TObjArray* fEMCALRecalibrationFactors; // Array of histograms with map of recalibration factors, EMCAL
361     
362   // Time Recalibration 
363   Bool_t     fTimeRecalibration;             // Switch on or off the time recalibration
364   TObjArray* fEMCALTimeRecalibrationFactors; // Array of histograms with map of time recalibration factors, EMCAL
365   
366   // Recalibrate with run dependent corrections, energy
367   Bool_t     fUseRunCorrectionFactors;   // Use Run Dependent Correction
368   Bool_t     fRunCorrectionFactorsSet;   // Run Correction set at leat once
369     
370   // Bad Channels
371   Bool_t     fRemoveBadChannels;         // Check the channel status provided and remove clusters with bad channels
372   Bool_t     fRecalDistToBadChannels;    // Calculate distance from highest energy tower of cluster to closes bad channel
373   TObjArray* fEMCALBadChannelMap;        // Array of histograms with map of bad channels, EMCAL
374
375   // Border cells
376   Int_t      fNCellsFromEMCALBorder;     // Number of cells from EMCAL border the cell with maximum amplitude has to be.
377   Bool_t     fNoEMCALBorderAtEta0;       // Do fiducial cut in EMCAL region eta = 0?
378   
379   // Exotic cell / cluster
380   Bool_t     fRejectExoticCluster;       // Switch on or off exotic cluster rejection
381   Bool_t     fRejectExoticCells;         // Remove exotic cells
382   Float_t    fExoticCellFraction;        // Good cell if fraction < 1-ecross/ecell
383   Float_t    fExoticCellDiffTime;        // If time of candidate to exotic and close cell is too different (in ns), it must be noisy, set amp to 0
384   Float_t    fExoticCellMinAmplitude;    // Check for exotic only if amplitud is larger than this value
385   
386   // PID
387   AliEMCALPIDUtils * fPIDUtils;          // Recalculate PID parameters
388     
389   //Track matching
390   UInt_t     fAODFilterMask;             // Filter mask to select AOD tracks. Refer to $ALICE_ROOT/ANALYSIS/macros/AddTaskESDFilter.C
391   TArrayI  * fMatchedTrackIndex;         // Array that stores indexes of matched tracks      
392   TArrayI  * fMatchedClusterIndex;       // Array that stores indexes of matched clusters
393   TArrayF  * fResidualEta;               // Array that stores the residual eta
394   TArrayF  * fResidualPhi;               // Array that stores the residual phi
395   Bool_t     fCutEtaPhiSum;              // Place cut on sqrt(dEta^2+dPhi^2)
396   Bool_t     fCutEtaPhiSeparate;         // Cut on dEta and dPhi separately
397   Float_t    fCutR;                      // sqrt(dEta^2+dPhi^2) cut on matching
398   Float_t    fCutEta;                    // dEta cut on matching
399   Float_t    fCutPhi;                    // dPhi cut on matching
400   Double_t   fClusterWindow;             // Select clusters in the window to be matched
401   Double_t   fMass;                      // Mass hypothesis of the track
402   Double_t   fStepSurface;               // Length of step to extrapolate tracks to EMCal surface
403   Double_t   fStepCluster;               // Length of step to extrapolate tracks to clusters
404
405   // Track cuts  
406   Int_t      fTrackCutsType;             // Esd track cuts type for matching
407   Double_t   fCutMinTrackPt;             // Cut on track pT
408   Int_t      fCutMinNClusterTPC;         // Min number of tpc clusters
409   Int_t      fCutMinNClusterITS;         // Min number of its clusters  
410   Float_t    fCutMaxChi2PerClusterTPC;   // Max tpc fit chi2 per tpc cluster
411   Float_t    fCutMaxChi2PerClusterITS;   // Max its fit chi2 per its cluster
412   Bool_t     fCutRequireTPCRefit;        // Require TPC refit
413   Bool_t     fCutRequireITSRefit;        // Require ITS refit
414   Bool_t     fCutAcceptKinkDaughters;    // Accepting kink daughters?
415   Float_t    fCutMaxDCAToVertexXY;       // Track-to-vertex cut in max absolute distance in xy-plane
416   Float_t    fCutMaxDCAToVertexZ;        // Track-to-vertex cut in max absolute distance in z-plane
417   Bool_t     fCutDCAToVertex2D;          // If true a 2D DCA cut is made.
418   
419   ClassDef(AliEMCALRecoUtils, 17)
420   
421 };
422
423 #endif // ALIEMCALRECOUTILS_H
424
425