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