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