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