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 ///////////////////////////////////////////////////////////////////////////////
34 class AliEMCALGeometry;
35 class AliEMCALPIDUtils;
37 class AliExternalTrackParam;
39 class AliEMCALRecoUtils : public TNamed {
44 AliEMCALRecoUtils(const AliEMCALRecoUtils&);
45 AliEMCALRecoUtils& operator=(const AliEMCALRecoUtils&);
46 virtual ~AliEMCALRecoUtils() ;
48 void InitParameters();
50 void Print(const Option_t*) const;
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};
59 //-----------------------------------------------------
60 //Position recalculation
61 //-----------------------------------------------------
63 void RecalculateClusterPosition (AliEMCALGeometry *geom, AliVCaloCells* cells, AliVCluster* clu);
64 void RecalculateClusterPositionFromTowerIndex (AliEMCALGeometry *geom, AliVCaloCells* cells, AliVCluster* clu);
65 void RecalculateClusterPositionFromTowerGlobal(AliEMCALGeometry *geom, AliVCaloCells* cells, AliVCluster* clu);
67 Float_t GetCellWeight(const Float_t eCell, const Float_t eCluster) const { return TMath::Max( 0., fW0 + TMath::Log( eCell / eCluster )) ; }
69 Float_t GetDepth(const Float_t eCluster, const Int_t iParticle, const Int_t iSM) const ;
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);
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)) ;
77 Float_t* GetMisalTransShiftArray() { return fMisalTransShift ; }
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] ; }
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)) ;
88 Float_t* GetMisalRotShiftArray() { return fMisalRotShift ; }
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)) ; } }
94 void SetMisalRotShiftArray(Float_t * misal) { for(Int_t i = 0; i < 15; i++)fMisalRotShift[i] = misal[i] ; }
96 Int_t GetParticleType() const { return fParticleType ; }
97 void SetParticleType(Int_t particle) { fParticleType = particle ; }
99 Int_t GetPositionAlgorithm() const { return fPosAlgo ; }
100 void SetPositionAlgorithm(Int_t alg) { fPosAlgo = alg ; }
102 Float_t GetW0() const { return fW0 ; }
103 void SetW0(Float_t w0) { fW0 = w0 ; }
105 //-----------------------------------------------------
107 //-----------------------------------------------------
109 Float_t CorrectClusterEnergyLinearity(AliVCluster* clu) ;
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)) ;
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();
119 Int_t GetNonLinearityFunction() const { return fNonLinearityFunction ; }
120 void SetNonLinearityFunction(Int_t fun) { fNonLinearityFunction = fun ; InitNonLinearityParam() ; }
122 void SetNonLinearityThreshold(Int_t threshold) { fNonLinearThreshold = threshold ; } //only for Alexie's non linearity correction
123 Int_t GetNonLinearityThreshold() const { return fNonLinearThreshold ; }
125 //-----------------------------------------------------
126 // MC clusters energy smearing
127 //-----------------------------------------------------
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 //-----------------------------------------------------
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
143 // Energy recalibration
144 Bool_t IsRecalibrationOn() const { return fRecalibration ; }
145 void SwitchOffRecalibration() { fRecalibration = kFALSE ; }
146 void SwitchOnRecalibration() { fRecalibration = kTRUE ;
147 if(!fEMCALRecalibrationFactors)InitEMCALRecalibrationFactors() ; }
148 void InitEMCALRecalibrationFactors() ;
150 TH2F * GetEMCALChannelRecalibrationFactors(Int_t iSM) const { return (TH2F*)fEMCALRecalibrationFactors->At(iSM) ; }
151 void SetEMCALChannelRecalibrationFactors(TObjArray *map) { fEMCALRecalibrationFactors = map ; }
152 void SetEMCALChannelRecalibrationFactors(Int_t iSM , TH2F* h) { fEMCALRecalibrationFactors->AddAt(h,iSM) ; }
154 Float_t GetEMCALChannelRecalibrationFactor(Int_t iSM , Int_t iCol, Int_t iRow) const {
155 if(fEMCALRecalibrationFactors)
156 return (Float_t) ((TH2F*)fEMCALRecalibrationFactors->At(iSM))->GetBinContent(iCol,iRow);
159 void SetEMCALChannelRecalibrationFactor(Int_t iSM , Int_t iCol, Int_t iRow, Double_t c = 1) {
160 if(!fEMCALRecalibrationFactors) InitEMCALRecalibrationFactors() ;
161 ((TH2F*)fEMCALRecalibrationFactors->At(iSM))->SetBinContent(iCol,iRow,c) ; }
163 //Recalibrate channels energy with run dependent corrections
164 void SwitchOffRunDepCorrection() { fUseRunCorrectionFactors = kFALSE ; }
165 void SwitchOnRunDepCorrection() { fUseRunCorrectionFactors = kTRUE ;
166 SwitchOnRecalibration() ; }
167 void SetRunDependentCorrections(Int_t runnumber);
169 // Time Recalibration
170 void RecalibrateCellTime(const Int_t absId, const Int_t bc, Double_t & time) const;
172 Bool_t IsTimeRecalibrationOn() const { return fTimeRecalibration ; }
173 void SwitchOffTimeRecalibration() { fTimeRecalibration = kFALSE ; }
174 void SwitchOnTimeRecalibration() { fTimeRecalibration = kTRUE ;
175 if(!fEMCALTimeRecalibrationFactors)InitEMCALTimeRecalibrationFactors() ; }
176 void InitEMCALTimeRecalibrationFactors() ;
178 Float_t GetEMCALChannelTimeRecalibrationFactor(const Int_t bc, const Int_t absID) const {
179 if(fEMCALTimeRecalibrationFactors)
180 return (Float_t) ((TH1F*)fEMCALTimeRecalibrationFactors->At(bc))->GetBinContent(absID);
183 void SetEMCALChannelTimeRecalibrationFactor(const Int_t bc, const Int_t absID, Double_t c = 0) {
184 if(!fEMCALTimeRecalibrationFactors) InitEMCALTimeRecalibrationFactors() ;
185 ((TH1F*)fEMCALTimeRecalibrationFactors->At(bc))->SetBinContent(absID,c) ; }
187 TH1F * GetEMCALChannelTimeRecalibrationFactors(const Int_t bc)const { return (TH1F*)fEMCALTimeRecalibrationFactors->At(bc) ; }
188 void SetEMCALChannelTimeRecalibrationFactors(TObjArray *map) { fEMCALTimeRecalibrationFactors = map ; }
189 void SetEMCALChannelTimeRecalibrationFactors(const Int_t bc , TH1F* h) { fEMCALTimeRecalibrationFactors->AddAt(h,bc) ; }
191 //-----------------------------------------------------
192 // Modules fiducial region, remove clusters in borders
193 //-----------------------------------------------------
195 Bool_t CheckCellFiducialRegion(AliEMCALGeometry* geom, AliVCluster* cluster, AliVCaloCells* cells) ;
196 void SetNumberOfCellsFromEMCALBorder(const Int_t n){ fNCellsFromEMCALBorder = n ; }
197 Int_t GetNumberOfCellsFromEMCALBorder() const { return fNCellsFromEMCALBorder ; }
199 void SwitchOnNoFiducialBorderInEMCALEta0() { fNoEMCALBorderAtEta0 = kTRUE ; }
200 void SwitchOffNoFiducialBorderInEMCALEta0() { fNoEMCALBorderAtEta0 = kFALSE ; }
201 Bool_t IsEMCALNoBorderAtEta0() const { return fNoEMCALBorderAtEta0 ; }
203 //-----------------------------------------------------
205 //-----------------------------------------------------
207 Bool_t IsBadChannelsRemovalSwitchedOn() const { return fRemoveBadChannels ; }
208 void SwitchOffBadChannelsRemoval() { fRemoveBadChannels = kFALSE ; }
209 void SwitchOnBadChannelsRemoval () { fRemoveBadChannels = kTRUE ;
210 if(!fEMCALBadChannelMap)InitEMCALBadChannelStatusMap() ; }
212 Bool_t IsDistanceToBadChannelRecalculated() const { return fRecalDistToBadChannels ; }
213 void SwitchOffDistToBadChannelRecalculation() { fRecalDistToBadChannels = kFALSE ; }
214 void SwitchOnDistToBadChannelRecalculation() { fRecalDistToBadChannels = kTRUE ;
215 if(!fEMCALBadChannelMap)InitEMCALBadChannelStatusMap() ; }
217 void InitEMCALBadChannelStatusMap() ;
219 Int_t GetEMCALChannelStatus(Int_t iSM , Int_t iCol, Int_t iRow) const {
220 if(fEMCALBadChannelMap) return (Int_t) ((TH2I*)fEMCALBadChannelMap->At(iSM))->GetBinContent(iCol,iRow);
221 else return 0;}//Channel is ok by default
223 void SetEMCALChannelStatus(Int_t iSM , Int_t iCol, Int_t iRow, Double_t c = 1) {
224 if(!fEMCALBadChannelMap)InitEMCALBadChannelStatusMap() ;
225 ((TH2I*)fEMCALBadChannelMap->At(iSM))->SetBinContent(iCol,iRow,c) ; }
227 TH2I * GetEMCALChannelStatusMap(Int_t iSM) const { return (TH2I*)fEMCALBadChannelMap->At(iSM) ; }
228 void SetEMCALChannelStatusMap(TObjArray *map) { fEMCALBadChannelMap = map ; }
229 void SetEMCALChannelStatusMap(Int_t iSM , TH2I* h) { fEMCALBadChannelMap->AddAt(h,iSM) ; }
231 Bool_t ClusterContainsBadChannel(const AliEMCALGeometry* geom, const UShort_t* cellList, const Int_t nCells);
233 //-----------------------------------------------------
234 // Recalculate other cluster parameters
235 //-----------------------------------------------------
237 void RecalculateClusterDistanceToBadChannel (AliEMCALGeometry * geom, AliVCaloCells* cells, AliVCluster * cluster);
238 void RecalculateClusterShowerShapeParameters(AliEMCALGeometry * geom, AliVCaloCells* cells, AliVCluster * cluster);
239 void RecalculateClusterPID(AliVCluster * cluster);
241 AliEMCALPIDUtils * GetPIDUtils() { return fPIDUtils;}
244 //----------------------------------------------------
246 //----------------------------------------------------
248 void FindMatches(AliVEvent *event, TObjArray * clusterArr=0x0, AliEMCALGeometry *geom=0x0);
249 Int_t FindMatchedClusterInEvent(AliESDtrack *track, AliVEvent *event, AliEMCALGeometry *geom, Float_t &dEta, Float_t &dPhi);
250 Int_t FindMatchedClusterInClusterArr(AliExternalTrackParam *emcalParam, AliExternalTrackParam *trkParam, TObjArray * clusterArr, Float_t &dEta, Float_t &dPhi);
252 static Bool_t ExtrapolateTrackToEMCalSurface(AliExternalTrackParam *trkParam, Double_t emcalR,
253 Double_t mass, Double_t step, Float_t &eta, Float_t &phi);
254 static Bool_t ExtrapolateTrackToPosition(AliExternalTrackParam *trkParam, const Float_t *clsPos,
255 Double_t mass, Double_t step, Float_t &tmpEta, Float_t &tmpPhi);
256 static Bool_t ExtrapolateTrackToCluster (AliExternalTrackParam *trkParam, AliVCluster *cluster,
257 Double_t mass, Double_t step, Float_t &tmpEta, Float_t &tmpPhi);
258 Bool_t ExtrapolateTrackToCluster (AliExternalTrackParam *trkParam, AliVCluster *cluster,
259 Float_t &tmpEta, Float_t &tmpPhi);
261 UInt_t FindMatchedPosForCluster(Int_t clsIndex) const;
262 UInt_t FindMatchedPosForTrack(Int_t trkIndex) const;
264 void GetMatchedResiduals(Int_t clsIndex, Float_t &dEta, Float_t &dPhi);
265 void GetMatchedClusterResiduals(Int_t trkIndex, Float_t &dEta, Float_t &dPhi);
266 Int_t GetMatchedTrackIndex(Int_t clsIndex);
267 Int_t GetMatchedClusterIndex(Int_t trkIndex);
269 Bool_t IsClusterMatched(Int_t clsIndex) const;
270 Bool_t IsTrackMatched(Int_t trkIndex) const;
272 void SetClusterMatchedToTrack (const AliESDEvent *event);
274 void SetTracksMatchedToCluster(const AliESDEvent *event);
276 void SwitchOnCutEtaPhiSum() { fCutEtaPhiSum = kTRUE ;
277 fCutEtaPhiSeparate = kFALSE ; }
278 void SwitchOnCutEtaPhiSeparate() { fCutEtaPhiSeparate = kTRUE ;
279 fCutEtaPhiSum = kFALSE ; }
281 Float_t GetCutR() const { return fCutR ; }
282 Float_t GetCutEta() const { return fCutEta ; }
283 Float_t GetCutPhi() const { return fCutPhi ; }
284 Double_t GetClusterWindow() const { return fClusterWindow ; }
285 void SetCutR(Float_t cutR) { fCutR = cutR ; }
286 void SetCutEta(Float_t cutEta) { fCutEta = cutEta ; }
287 void SetCutPhi(Float_t cutPhi) { fCutPhi = cutPhi ; }
288 void SetClusterWindow(Double_t window) { fClusterWindow = window ; }
289 void SetCutZ(Float_t cutZ) { printf("Obsolete fucntion of cutZ=%1.1f\n",cutZ) ; } //Obsolete
291 Double_t GetMass() const { return fMass ; }
292 Double_t GetStep() const { return fStepCluster ; }
293 Double_t GetStepSurface() const { return fStepSurface ; }
294 void SetMass(Double_t mass) { fMass = mass ; }
295 void SetStep(Double_t step) { fStepSurface = step ; }
296 void SetStepCluster(Double_t step) { fStepCluster = step ; }
298 // Exotic cells / clusters
300 Bool_t IsExoticCell(const Int_t absId, AliVCaloCells* cells, const Int_t bc =-1) ;
301 void SwitchOnRejectExoticCell() { fRejectExoticCells = kTRUE ; }
302 void SwitchOffRejectExoticCell() { fRejectExoticCells = kFALSE ; }
304 void SetExoticCellFractionCut(Float_t f) { fExoticCellFraction = f ; }
305 void SetExoticCellDiffTimeCut(Float_t dt) { fExoticCellDiffTime = dt ; }
306 void SetExoticCellMinAmplitudeCut(Float_t ma) { fExoticCellMinAmplitude = ma ; }
308 Bool_t IsExoticCluster(AliVCluster *cluster, AliVCaloCells* cells, const Int_t bc=0) ;
309 void SwitchOnRejectExoticCluster() { fRejectExoticCluster = kTRUE ;
310 fRejectExoticCells = kTRUE ; }
311 void SwitchOffRejectExoticCluster() { fRejectExoticCluster = kFALSE ; }
312 Bool_t IsRejectExoticCluster() const { return fRejectExoticCluster ; }
315 Bool_t IsGoodCluster(AliVCluster *cluster, AliEMCALGeometry *geom, AliVCaloCells* cells, const Int_t bc =-1);
318 Bool_t IsAccepted(AliESDtrack *track);
319 void InitTrackCuts();
320 void SetTrackCutsType(Int_t type) { fTrackCutsType = type ;
322 Int_t GetTrackCutsType() const { return fTrackCutsType; }
324 // track quality cut setters
325 void SetMinTrackPt(Double_t pt=0) { fCutMinTrackPt = pt ; }
326 void SetMinNClustersTPC(Int_t min=-1) { fCutMinNClusterTPC = min ; }
327 void SetMinNClustersITS(Int_t min=-1) { fCutMinNClusterITS = min ; }
328 void SetMaxChi2PerClusterTPC(Float_t max=1e10) { fCutMaxChi2PerClusterTPC = max ; }
329 void SetMaxChi2PerClusterITS(Float_t max=1e10) { fCutMaxChi2PerClusterITS = max ; }
330 void SetRequireTPCRefit(Bool_t b=kFALSE) { fCutRequireTPCRefit = b ; }
331 void SetRequireITSRefit(Bool_t b=kFALSE) { fCutRequireITSRefit = b ; }
332 void SetAcceptKinkDaughters(Bool_t b=kTRUE) { fCutAcceptKinkDaughters = b ; }
333 void SetMaxDCAToVertexXY(Float_t dist=1e10) { fCutMaxDCAToVertexXY = dist ; }
334 void SetMaxDCAToVertexZ(Float_t dist=1e10) { fCutMaxDCAToVertexZ = dist ; }
335 void SetDCAToVertex2D(Bool_t b=kFALSE) { fCutDCAToVertex2D = b ; }
338 Double_t GetMinTrackPt() const { return fCutMinTrackPt ; }
339 Int_t GetMinNClusterTPC() const { return fCutMinNClusterTPC ; }
340 Int_t GetMinNClustersITS() const { return fCutMinNClusterITS ; }
341 Float_t GetMaxChi2PerClusterTPC() const { return fCutMaxChi2PerClusterTPC ; }
342 Float_t GetMaxChi2PerClusterITS() const { return fCutMaxChi2PerClusterITS ; }
343 Bool_t GetRequireTPCRefit() const { return fCutRequireTPCRefit ; }
344 Bool_t GetRequireITSRefit() const { return fCutRequireITSRefit ; }
345 Bool_t GetAcceptKinkDaughters() const { return fCutAcceptKinkDaughters ; }
346 Float_t GetMaxDCAToVertexXY() const { return fCutMaxDCAToVertexXY ; }
347 Float_t GetMaxDCAToVertexZ() const { return fCutMaxDCAToVertexZ ; }
348 Bool_t GetDCAToVertex2D() const { return fCutDCAToVertex2D ; }
352 //Position recalculation
353 Float_t fMisalTransShift[15]; // Shift parameters
354 Float_t fMisalRotShift[15]; // Shift parameters
355 Int_t fParticleType; // Particle type for depth calculation
356 Int_t fPosAlgo; // Position recalculation algorithm
357 Float_t fW0; // Weight0
360 Int_t fNonLinearityFunction; // Non linearity function choice
361 Float_t fNonLinearityParams[7]; // Parameters for the non linearity function
362 Int_t fNonLinearThreshold; // Non linearity threshold value for kBeamTesh non linearity function
364 // Energy smearing for MC
365 Bool_t fSmearClusterEnergy; // Smear cluster energy, to be done only for simulated data to match real data
366 Float_t fSmearClusterParam[3]; // Smearing parameters
367 TRandom3 fRandom; // Random generator
369 // Energy Recalibration
370 Bool_t fCellsRecalibrated; // Internal bool to check if cells (time/energy) where recalibrated and not recalibrate them when recalculating different things
371 Bool_t fRecalibration; // Switch on or off the recalibration
372 TObjArray* fEMCALRecalibrationFactors; // Array of histograms with map of recalibration factors, EMCAL
374 // Time Recalibration
375 Bool_t fTimeRecalibration; // Switch on or off the time recalibration
376 TObjArray* fEMCALTimeRecalibrationFactors; // Array of histograms with map of time recalibration factors, EMCAL
378 // Recalibrate with run dependent corrections, energy
379 Bool_t fUseRunCorrectionFactors; // Use Run Dependent Correction
380 Bool_t fRunCorrectionFactorsSet; // Run Correction set at leat once
383 Bool_t fRemoveBadChannels; // Check the channel status provided and remove clusters with bad channels
384 Bool_t fRecalDistToBadChannels; // Calculate distance from highest energy tower of cluster to closes bad channel
385 TObjArray* fEMCALBadChannelMap; // Array of histograms with map of bad channels, EMCAL
388 Int_t fNCellsFromEMCALBorder; // Number of cells from EMCAL border the cell with maximum amplitude has to be.
389 Bool_t fNoEMCALBorderAtEta0; // Do fiducial cut in EMCAL region eta = 0?
391 // Exotic cell / cluster
392 Bool_t fRejectExoticCluster; // Switch on or off exotic cluster rejection
393 Bool_t fRejectExoticCells; // Remove exotic cells
394 Float_t fExoticCellFraction; // Good cell if fraction < 1-ecross/ecell
395 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
396 Float_t fExoticCellMinAmplitude; // Check for exotic only if amplitud is larger than this value
399 AliEMCALPIDUtils * fPIDUtils; // Recalculate PID parameters
402 UInt_t fAODFilterMask; // Filter mask to select AOD tracks. Refer to $ALICE_ROOT/ANALYSIS/macros/AddTaskESDFilter.C
403 TArrayI * fMatchedTrackIndex; // Array that stores indexes of matched tracks
404 TArrayI * fMatchedClusterIndex; // Array that stores indexes of matched clusters
405 TArrayF * fResidualEta; // Array that stores the residual eta
406 TArrayF * fResidualPhi; // Array that stores the residual phi
407 Bool_t fCutEtaPhiSum; // Place cut on sqrt(dEta^2+dPhi^2)
408 Bool_t fCutEtaPhiSeparate; // Cut on dEta and dPhi separately
409 Float_t fCutR; // sqrt(dEta^2+dPhi^2) cut on matching
410 Float_t fCutEta; // dEta cut on matching
411 Float_t fCutPhi; // dPhi cut on matching
412 Double_t fClusterWindow; // Select clusters in the window to be matched
413 Double_t fMass; // Mass hypothesis of the track
414 Double_t fStepSurface; // Length of step to extrapolate tracks to EMCal surface
415 Double_t fStepCluster; // Length of step to extrapolate tracks to clusters
418 Int_t fTrackCutsType; // Esd track cuts type for matching
419 Double_t fCutMinTrackPt; // Cut on track pT
420 Int_t fCutMinNClusterTPC; // Min number of tpc clusters
421 Int_t fCutMinNClusterITS; // Min number of its clusters
422 Float_t fCutMaxChi2PerClusterTPC; // Max tpc fit chi2 per tpc cluster
423 Float_t fCutMaxChi2PerClusterITS; // Max its fit chi2 per its cluster
424 Bool_t fCutRequireTPCRefit; // Require TPC refit
425 Bool_t fCutRequireITSRefit; // Require ITS refit
426 Bool_t fCutAcceptKinkDaughters; // Accepting kink daughters?
427 Float_t fCutMaxDCAToVertexXY; // Track-to-vertex cut in max absolute distance in xy-plane
428 Float_t fCutMaxDCAToVertexZ; // Track-to-vertex cut in max absolute distance in z-plane
429 Bool_t fCutDCAToVertex2D; // If true a 2D DCA cut is made.
431 ClassDef(AliEMCALRecoUtils, 17)
435 #endif // ALIEMCALRECOUTILS_H