1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
16 /* $Id: AliEMCALRecoUtils.cxx 33808 2009-07-15 09:48:08Z gconesab $ */
18 ///////////////////////////////////////////////////////////////////////////////
20 // Class AliEMCALRecoUtils
21 // Some utilities to recalculate the cluster position or energy linearity
24 // Author: Gustavo Conesa (LPSC- Grenoble)
25 // Track matching part: Rongrong Ma (Yale)
27 ///////////////////////////////////////////////////////////////////////////////
30 // standard C++ includes
31 //#include <Riostream.h>
34 #include <TGeoManager.h>
35 #include <TGeoMatrix.h>
40 #include <TObjArray.h>
43 #include "AliVCluster.h"
44 #include "AliVCaloCells.h"
47 #include "AliESDEvent.h"
48 #include "AliAODEvent.h"
49 #include "AliESDtrack.h"
50 #include "AliAODTrack.h"
51 #include "AliExternalTrackParam.h"
52 #include "AliESDfriendTrack.h"
53 #include "AliTrackerBase.h"
56 #include "AliEMCALRecoUtils.h"
57 #include "AliEMCALGeometry.h"
58 #include "AliTrackerBase.h"
59 #include "AliEMCALPIDUtils.h"
61 #include "AliESDtrackCuts.h" //MARCEL
62 ClassImp(AliEMCALRecoUtils)
64 //_____________________________________
65 AliEMCALRecoUtils::AliEMCALRecoUtils():
66 fParticleType(0), fPosAlgo(0), fW0(0),
67 fNonLinearityFunction(0), fNonLinearThreshold(0),
68 fSmearClusterEnergy(kFALSE), fRandom(),
69 fCellsRecalibrated(kFALSE), fRecalibration(kFALSE), fEMCALRecalibrationFactors(),
70 fTimeRecalibration(kFALSE), fEMCALTimeRecalibrationFactors(), fUseRunCorrectionFactors(kFALSE),
71 fRemoveBadChannels(kFALSE), fRecalDistToBadChannels(kFALSE), fEMCALBadChannelMap(),
72 fNCellsFromEMCALBorder(0), fNoEMCALBorderAtEta0(kTRUE),
73 fRejectExoticCluster(kFALSE), fRejectExoticCells(kFALSE),
74 fExoticCellFraction(0), fExoticCellDiffTime(0), fExoticCellMinAmplitude(0),
75 fPIDUtils(), fAODFilterMask(0),
76 fMatchedTrackIndex(0x0), fMatchedClusterIndex(0x0),
77 fResidualEta(0x0), fResidualPhi(0x0), fCutEtaPhiSum(kFALSE), fCutEtaPhiSeparate(kFALSE),
78 fCutR(0), fCutEta(0), fCutPhi(0),
79 fClusterWindow(0), fMass(0),
80 fStepSurface(0), fStepCluster(0),
82 fTrackCutsType(0), fCutMinTrackPt(0), fCutMinNClusterTPC(0),
83 fCutMinNClusterITS(0), fCutMaxChi2PerClusterTPC(0), fCutMaxChi2PerClusterITS(0),
84 fCutRequireTPCRefit(kFALSE), fCutRequireITSRefit(kFALSE), fCutAcceptKinkDaughters(kFALSE),
85 fCutMaxDCAToVertexXY(0), fCutMaxDCAToVertexZ(0), fCutDCAToVertex2D(kFALSE),
86 fCutRequireITSStandAlone(kFALSE), fCutRequireITSpureSA(kFALSE)
90 // Initialize all constant values which have to be used
91 // during Reco algorithm execution
98 fMatchedTrackIndex = new TArrayI();
99 fMatchedClusterIndex = new TArrayI();
100 fResidualPhi = new TArrayF();
101 fResidualEta = new TArrayF();
102 fPIDUtils = new AliEMCALPIDUtils();
106 //______________________________________________________________________
107 AliEMCALRecoUtils::AliEMCALRecoUtils(const AliEMCALRecoUtils & reco)
109 fParticleType(reco.fParticleType), fPosAlgo(reco.fPosAlgo), fW0(reco.fW0),
110 fNonLinearityFunction(reco.fNonLinearityFunction), fNonLinearThreshold(reco.fNonLinearThreshold),
111 fSmearClusterEnergy(reco.fSmearClusterEnergy), fRandom(),
112 fCellsRecalibrated(reco.fCellsRecalibrated),
113 fRecalibration(reco.fRecalibration), fEMCALRecalibrationFactors(reco.fEMCALRecalibrationFactors),
114 fTimeRecalibration(reco.fTimeRecalibration), fEMCALTimeRecalibrationFactors(reco.fEMCALTimeRecalibrationFactors),
115 fUseRunCorrectionFactors(reco.fUseRunCorrectionFactors),
116 fRemoveBadChannels(reco.fRemoveBadChannels), fRecalDistToBadChannels(reco.fRecalDistToBadChannels),
117 fEMCALBadChannelMap(reco.fEMCALBadChannelMap),
118 fNCellsFromEMCALBorder(reco.fNCellsFromEMCALBorder), fNoEMCALBorderAtEta0(reco.fNoEMCALBorderAtEta0),
119 fRejectExoticCluster(reco.fRejectExoticCluster), fRejectExoticCells(reco.fRejectExoticCells),
120 fExoticCellFraction(reco.fExoticCellFraction), fExoticCellDiffTime(reco.fExoticCellDiffTime),
121 fExoticCellMinAmplitude(reco.fExoticCellMinAmplitude),
122 fPIDUtils(reco.fPIDUtils), fAODFilterMask(reco.fAODFilterMask),
123 fMatchedTrackIndex( reco.fMatchedTrackIndex? new TArrayI(*reco.fMatchedTrackIndex):0x0),
124 fMatchedClusterIndex(reco.fMatchedClusterIndex?new TArrayI(*reco.fMatchedClusterIndex):0x0),
125 fResidualEta( reco.fResidualEta? new TArrayF(*reco.fResidualEta):0x0),
126 fResidualPhi( reco.fResidualPhi? new TArrayF(*reco.fResidualPhi):0x0),
127 fCutEtaPhiSum(reco.fCutEtaPhiSum), fCutEtaPhiSeparate(reco.fCutEtaPhiSeparate),
128 fCutR(reco.fCutR), fCutEta(reco.fCutEta), fCutPhi(reco.fCutPhi),
129 fClusterWindow(reco.fClusterWindow),
130 fMass(reco.fMass), fStepSurface(reco.fStepSurface), fStepCluster(reco.fStepCluster),
131 fITSTrackSA(reco.fITSTrackSA),
132 fTrackCutsType(reco.fTrackCutsType), fCutMinTrackPt(reco.fCutMinTrackPt),
133 fCutMinNClusterTPC(reco.fCutMinNClusterTPC), fCutMinNClusterITS(reco.fCutMinNClusterITS),
134 fCutMaxChi2PerClusterTPC(reco.fCutMaxChi2PerClusterTPC), fCutMaxChi2PerClusterITS(reco.fCutMaxChi2PerClusterITS),
135 fCutRequireTPCRefit(reco.fCutRequireTPCRefit), fCutRequireITSRefit(reco.fCutRequireITSRefit),
136 fCutAcceptKinkDaughters(reco.fCutAcceptKinkDaughters), fCutMaxDCAToVertexXY(reco.fCutMaxDCAToVertexXY),
137 fCutMaxDCAToVertexZ(reco.fCutMaxDCAToVertexZ), fCutDCAToVertex2D(reco.fCutDCAToVertex2D),
138 fCutRequireITSStandAlone(reco.fCutRequireITSStandAlone), fCutRequireITSpureSA(reco.fCutRequireITSpureSA)
142 for(Int_t i = 0; i < 15 ; i++) { fMisalRotShift[i] = reco.fMisalRotShift[i] ;
143 fMisalTransShift[i] = reco.fMisalTransShift[i] ; }
144 for(Int_t i = 0; i < 7 ; i++) { fNonLinearityParams[i] = reco.fNonLinearityParams[i] ; }
145 for(Int_t i = 0; i < 3 ; i++) { fSmearClusterParam[i] = reco.fSmearClusterParam[i] ; }
150 //______________________________________________________________________
151 AliEMCALRecoUtils & AliEMCALRecoUtils::operator = (const AliEMCALRecoUtils & reco)
153 //Assignment operator
155 if(this == &reco)return *this;
156 ((TNamed *)this)->operator=(reco);
158 for(Int_t i = 0; i < 15 ; i++) { fMisalTransShift[i] = reco.fMisalTransShift[i] ;
159 fMisalRotShift[i] = reco.fMisalRotShift[i] ; }
160 for(Int_t i = 0; i < 7 ; i++) { fNonLinearityParams[i] = reco.fNonLinearityParams[i] ; }
161 for(Int_t i = 0; i < 3 ; i++) { fSmearClusterParam[i] = reco.fSmearClusterParam[i] ; }
163 fParticleType = reco.fParticleType;
164 fPosAlgo = reco.fPosAlgo;
167 fNonLinearityFunction = reco.fNonLinearityFunction;
168 fNonLinearThreshold = reco.fNonLinearThreshold;
169 fSmearClusterEnergy = reco.fSmearClusterEnergy;
171 fCellsRecalibrated = reco.fCellsRecalibrated;
172 fRecalibration = reco.fRecalibration;
173 fEMCALRecalibrationFactors = reco.fEMCALRecalibrationFactors;
175 fTimeRecalibration = reco.fTimeRecalibration;
176 fEMCALTimeRecalibrationFactors = reco.fEMCALTimeRecalibrationFactors;
178 fUseRunCorrectionFactors = reco.fUseRunCorrectionFactors;
180 fRemoveBadChannels = reco.fRemoveBadChannels;
181 fRecalDistToBadChannels = reco.fRecalDistToBadChannels;
182 fEMCALBadChannelMap = reco.fEMCALBadChannelMap;
184 fNCellsFromEMCALBorder = reco.fNCellsFromEMCALBorder;
185 fNoEMCALBorderAtEta0 = reco.fNoEMCALBorderAtEta0;
187 fRejectExoticCluster = reco.fRejectExoticCluster;
188 fRejectExoticCells = reco.fRejectExoticCells;
189 fExoticCellFraction = reco.fExoticCellFraction;
190 fExoticCellDiffTime = reco.fExoticCellDiffTime;
191 fExoticCellMinAmplitude = reco.fExoticCellMinAmplitude;
193 fPIDUtils = reco.fPIDUtils;
195 fAODFilterMask = reco.fAODFilterMask;
197 fCutEtaPhiSum = reco.fCutEtaPhiSum;
198 fCutEtaPhiSeparate = reco.fCutEtaPhiSeparate;
200 fCutEta = reco.fCutEta;
201 fCutPhi = reco.fCutPhi;
202 fClusterWindow = reco.fClusterWindow;
204 fStepSurface = reco.fStepSurface;
205 fStepCluster = reco.fStepCluster;
206 fITSTrackSA = reco.fITSTrackSA;
208 fTrackCutsType = reco.fTrackCutsType;
209 fCutMinTrackPt = reco.fCutMinTrackPt;
210 fCutMinNClusterTPC = reco.fCutMinNClusterTPC;
211 fCutMinNClusterITS = reco.fCutMinNClusterITS;
212 fCutMaxChi2PerClusterTPC = reco.fCutMaxChi2PerClusterTPC;
213 fCutMaxChi2PerClusterITS = reco.fCutMaxChi2PerClusterITS;
214 fCutRequireTPCRefit = reco.fCutRequireTPCRefit;
215 fCutRequireITSRefit = reco.fCutRequireITSRefit;
216 fCutAcceptKinkDaughters = reco.fCutAcceptKinkDaughters;
217 fCutMaxDCAToVertexXY = reco.fCutMaxDCAToVertexXY;
218 fCutMaxDCAToVertexZ = reco.fCutMaxDCAToVertexZ;
219 fCutDCAToVertex2D = reco.fCutDCAToVertex2D;
220 fCutRequireITSStandAlone = reco.fCutRequireITSStandAlone;
221 fCutRequireITSpureSA = reco.fCutRequireITSpureSA;
222 if(reco.fResidualEta)
224 // assign or copy construct
227 *fResidualEta = *reco.fResidualEta;
231 fResidualEta = new TArrayF(*reco.fResidualEta);
236 if(fResidualEta)delete fResidualEta;
240 if(reco.fResidualPhi)
242 // assign or copy construct
245 *fResidualPhi = *reco.fResidualPhi;
249 fResidualPhi = new TArrayF(*reco.fResidualPhi);
254 if(fResidualPhi)delete fResidualPhi;
258 if(reco.fMatchedTrackIndex)
260 // assign or copy construct
261 if(fMatchedTrackIndex)
263 *fMatchedTrackIndex = *reco.fMatchedTrackIndex;
267 fMatchedTrackIndex = new TArrayI(*reco.fMatchedTrackIndex);
272 if(fMatchedTrackIndex)delete fMatchedTrackIndex;
273 fMatchedTrackIndex = 0;
276 if(reco.fMatchedClusterIndex)
278 // assign or copy construct
279 if(fMatchedClusterIndex)
281 *fMatchedClusterIndex = *reco.fMatchedClusterIndex;
285 fMatchedClusterIndex = new TArrayI(*reco.fMatchedClusterIndex);
290 if(fMatchedClusterIndex)delete fMatchedClusterIndex;
291 fMatchedClusterIndex = 0;
298 //_____________________________________
299 AliEMCALRecoUtils::~AliEMCALRecoUtils()
303 if(fEMCALRecalibrationFactors)
305 fEMCALRecalibrationFactors->Clear();
306 delete fEMCALRecalibrationFactors;
309 if(fEMCALTimeRecalibrationFactors)
311 fEMCALTimeRecalibrationFactors->Clear();
312 delete fEMCALTimeRecalibrationFactors;
315 if(fEMCALBadChannelMap)
317 fEMCALBadChannelMap->Clear();
318 delete fEMCALBadChannelMap;
321 delete fMatchedTrackIndex ;
322 delete fMatchedClusterIndex ;
323 delete fResidualEta ;
324 delete fResidualPhi ;
330 //_______________________________________________________________________________
331 Bool_t AliEMCALRecoUtils::AcceptCalibrateCell(const Int_t absID, const Int_t bc,
332 Float_t & amp, Double_t & time,
333 AliVCaloCells* cells)
335 // Reject cell if criteria not passed and calibrate it
337 AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
339 if(absID < 0 || absID >= 24*48*geom->GetNumberOfSuperModules()) return kFALSE;
341 Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1;
343 if(!geom->GetCellIndex(absID,imod,iTower,iIphi,iIeta))
345 // cell absID does not exist
350 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);
352 // Do not include bad channels found in analysis,
353 if( IsBadChannelsRemovalSwitchedOn() && GetEMCALChannelStatus(imod, ieta, iphi))
359 amp = cells->GetCellAmplitude(absID);
360 if(!fCellsRecalibrated && IsRecalibrationOn())
361 amp *= GetEMCALChannelRecalibrationFactor(imod,ieta,iphi);
365 time = cells->GetCellTime(absID);
367 RecalibrateCellTime(absID,bc,time);
372 //_____________________________________________________________________________
373 Bool_t AliEMCALRecoUtils::CheckCellFiducialRegion(const AliEMCALGeometry* geom,
374 const AliVCluster* cluster,
375 AliVCaloCells* cells)
377 // Given the list of AbsId of the cluster, get the maximum cell and
378 // check if there are fNCellsFromBorder from the calorimeter border
382 AliInfo("Cluster pointer null!");
386 //If the distance to the border is 0 or negative just exit accept all clusters
387 if(cells->GetType()==AliVCaloCells::kEMCALCell && fNCellsFromEMCALBorder <= 0 ) return kTRUE;
389 Int_t absIdMax = -1, iSM =-1, ieta = -1, iphi = -1;
390 Bool_t shared = kFALSE;
391 GetMaxEnergyCell(geom, cells, cluster, absIdMax, iSM, ieta, iphi, shared);
393 AliDebug(2,Form("Cluster Max AbsId %d, Cell Energy %2.2f, Cluster Energy %2.2f, Ncells from border %d, EMCAL eta=0 %d\n",
394 absIdMax, cells->GetCellAmplitude(absIdMax), cluster->E(), fNCellsFromEMCALBorder, fNoEMCALBorderAtEta0));
396 if(absIdMax==-1) return kFALSE;
398 //Check if the cell is close to the borders:
399 Bool_t okrow = kFALSE;
400 Bool_t okcol = kFALSE;
402 if(iSM < 0 || iphi < 0 || ieta < 0 )
404 AliFatal(Form("Negative value for super module: %d, or cell ieta: %d, or cell iphi: %d, check EMCAL geometry name\n",
411 if(iphi >= fNCellsFromEMCALBorder && iphi < 24-fNCellsFromEMCALBorder) okrow =kTRUE;
413 else if (iSM >=10 && ( ( geom->GetEMCGeometry()->GetGeoName()).Contains("12SMV1")))
415 if(iphi >= fNCellsFromEMCALBorder && iphi < 8-fNCellsFromEMCALBorder) okrow =kTRUE; //1/3 sm case
419 if(iphi >= fNCellsFromEMCALBorder && iphi < 12-fNCellsFromEMCALBorder) okrow =kTRUE; // half SM case
423 if(!fNoEMCALBorderAtEta0)
425 if(ieta > fNCellsFromEMCALBorder && ieta < 48-fNCellsFromEMCALBorder) okcol =kTRUE;
431 if(ieta >= fNCellsFromEMCALBorder) okcol = kTRUE;
435 if(ieta < 48-fNCellsFromEMCALBorder) okcol = kTRUE;
439 AliDebug(2,Form("EMCAL Cluster in %d cells fiducial volume: ieta %d, iphi %d, SM %d: column? %d, row? %d\nq",
440 fNCellsFromEMCALBorder, ieta, iphi, iSM, okcol, okrow));
444 //printf("Accept\n");
449 //printf("Reject\n");
450 AliDebug(2,Form("Reject cluster in border, max cell : ieta %d, iphi %d, SM %d\n",ieta, iphi, iSM));
457 //_______________________________________________________________________________
458 Bool_t AliEMCALRecoUtils::ClusterContainsBadChannel(const AliEMCALGeometry* geom,
459 const UShort_t* cellList,
462 // Check that in the cluster cells, there is no bad channel of those stored
463 // in fEMCALBadChannelMap or fPHOSBadChannelMap
465 if(!fRemoveBadChannels) return kFALSE;
466 if(!fEMCALBadChannelMap) return kFALSE;
471 for(Int_t iCell = 0; iCell<nCells; iCell++)
473 //Get the column and row
474 Int_t iTower = -1, iIphi = -1, iIeta = -1;
475 geom->GetCellIndex(cellList[iCell],imod,iTower,iIphi,iIeta);
476 if(fEMCALBadChannelMap->GetEntries() <= imod) continue;
477 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,irow,icol);
478 if(GetEMCALChannelStatus(imod, icol, irow))
480 AliDebug(2,Form("Cluster with bad channel: SM %d, col %d, row %d\n",imod, icol, irow));
484 }// cell cluster loop
490 //___________________________________________________________________________
491 Float_t AliEMCALRecoUtils::GetECross(const Int_t absID, const Double_t tcell,
492 AliVCaloCells* cells, const Int_t bc)
494 //Calculate the energy in the cross around the energy given cell
496 AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
498 Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1;
499 geom->GetCellIndex(absID,imod,iTower,iIphi,iIeta);
500 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);
502 //Get close cells index, energy and time, not in corners
507 if( iphi < AliEMCALGeoParams::fgkEMCALRows-1) absID1 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi+1, ieta);
508 if( iphi > 0 ) absID2 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi-1, ieta);
510 // In case of cell in eta = 0 border, depending on SM shift the cross cell index
515 if ( ieta == AliEMCALGeoParams::fgkEMCALCols-1 && !(imod%2) )
517 absID3 = geom-> GetAbsCellIdFromCellIndexes(imod+1, iphi, 0);
518 absID4 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta-1);
520 else if( ieta == 0 && imod%2 )
522 absID3 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta+1);
523 absID4 = geom-> GetAbsCellIdFromCellIndexes(imod-1, iphi, AliEMCALGeoParams::fgkEMCALCols-1);
527 if( ieta < AliEMCALGeoParams::fgkEMCALCols-1 )
528 absID3 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta+1);
530 absID4 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta-1);
533 //printf("IMOD %d, AbsId %d, a %d, b %d, c %d e %d \n",imod,absID,absID1,absID2,absID3,absID4);
535 Float_t ecell1 = 0, ecell2 = 0, ecell3 = 0, ecell4 = 0;
536 Double_t tcell1 = 0, tcell2 = 0, tcell3 = 0, tcell4 = 0;
537 Bool_t accept1 = 0, accept2 = 0, accept3 = 0, accept4 = 0;
539 accept1 = AcceptCalibrateCell(absID1,bc, ecell1,tcell1,cells);
540 accept2 = AcceptCalibrateCell(absID2,bc, ecell2,tcell2,cells);
541 accept3 = AcceptCalibrateCell(absID3,bc, ecell3,tcell3,cells);
542 accept4 = AcceptCalibrateCell(absID4,bc, ecell4,tcell4,cells);
545 printf("Cell absID %d \n",absID);
546 printf("\t accept1 %d, accept2 %d, accept3 %d, accept4 %d\n",
547 accept1,accept2,accept3,accept4);
548 printf("\t id %d: id1 %d, id2 %d, id3 %d, id4 %d\n",
549 absID,absID1,absID2,absID3,absID4);
550 printf("\t e %f: e1 %f, e2 %f, e3 %f, e4 %f\n",
551 ecell,ecell1,ecell2,ecell3,ecell4);
552 printf("\t t %f: t1 %f, t2 %f, t3 %f, t4 %f;\n dt1 %f, dt2 %f, dt3 %f, dt4 %f\n",
553 tcell*1.e9,tcell1*1.e9,tcell2*1.e9,tcell3*1.e9,tcell4*1.e9,
554 TMath::Abs(tcell-tcell1)*1.e9, TMath::Abs(tcell-tcell2)*1.e9, TMath::Abs(tcell-tcell3)*1.e9, TMath::Abs(tcell-tcell4)*1.e9);
557 if(TMath::Abs(tcell-tcell1)*1.e9 > fExoticCellDiffTime) ecell1 = 0 ;
558 if(TMath::Abs(tcell-tcell2)*1.e9 > fExoticCellDiffTime) ecell2 = 0 ;
559 if(TMath::Abs(tcell-tcell3)*1.e9 > fExoticCellDiffTime) ecell3 = 0 ;
560 if(TMath::Abs(tcell-tcell4)*1.e9 > fExoticCellDiffTime) ecell4 = 0 ;
562 return ecell1+ecell2+ecell3+ecell4;
566 //_____________________________________________________________________________________________
567 Bool_t AliEMCALRecoUtils::IsExoticCell(const Int_t absID, AliVCaloCells* cells, const Int_t bc)
569 // Look to cell neighbourhood and reject if it seems exotic
570 // Do before recalibrating the cells
572 if(!fRejectExoticCells) return kFALSE;
576 Bool_t accept = AcceptCalibrateCell(absID, bc, ecell ,tcell ,cells);
578 if(!accept) return kTRUE; // reject this cell
580 if(ecell < fExoticCellMinAmplitude) return kFALSE; // do not reject low energy cells
582 Float_t eCross = GetECross(absID,tcell,cells,bc);
584 if(1-eCross/ecell > fExoticCellFraction)
586 AliDebug(2,Form("AliEMCALRecoUtils::IsExoticCell() - EXOTIC CELL id %d, eCell %f, eCross %f, 1-eCross/eCell %f\n",
587 absID,ecell,eCross,1-eCross/ecell));
594 //___________________________________________________________________
595 Bool_t AliEMCALRecoUtils::IsExoticCluster(const AliVCluster *cluster,
596 AliVCaloCells *cells,
599 // Check if the cluster highest energy tower is exotic
603 AliInfo("Cluster pointer null!");
607 if(!fRejectExoticCluster) return kFALSE;
609 // Get highest energy tower
610 AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
611 Int_t iSupMod = -1, absId = -1, ieta = -1, iphi = -1;
612 Bool_t shared = kFALSE;
613 GetMaxEnergyCell(geom, cells, cluster, absId, iSupMod, ieta, iphi, shared);
615 return IsExoticCell(absId,cells,bc);
619 //_______________________________________________________________________
620 Float_t AliEMCALRecoUtils::SmearClusterEnergy(const AliVCluster* cluster)
622 //In case of MC analysis, smear energy to match resolution/calibration in real data
626 AliInfo("Cluster pointer null!");
630 Float_t energy = cluster->E() ;
631 Float_t rdmEnergy = energy ;
632 if(fSmearClusterEnergy)
634 rdmEnergy = fRandom.Gaus(energy,fSmearClusterParam[0] * TMath::Sqrt(energy) +
635 fSmearClusterParam[1] * energy +
636 fSmearClusterParam[2] );
637 AliDebug(2, Form("Energy: original %f, smeared %f\n", energy, rdmEnergy));
643 //____________________________________________________________________________
644 Float_t AliEMCALRecoUtils::CorrectClusterEnergyLinearity(AliVCluster* cluster)
646 // Correct cluster energy from non linearity functions
650 AliInfo("Cluster pointer null!");
654 Float_t energy = cluster->E();
658 // Clusters with less than 100 MeV or negative are not possible
659 AliInfo(Form("Too Low Cluster energy!, E = %f < 0.1 GeV",energy));
663 switch (fNonLinearityFunction)
668 //Non-Linearity correction (from MC with function ([0]*exp(-[1]/E))+(([2]/([3]*2.*TMath::Pi())*exp(-(E-[4])^2/(2.*[3]^2)))))
669 //fNonLinearityParams[0] = 1.014;
670 //fNonLinearityParams[1] =-0.03329;
671 //fNonLinearityParams[2] =-0.3853;
672 //fNonLinearityParams[3] = 0.5423;
673 //fNonLinearityParams[4] =-0.4335;
674 energy *= (fNonLinearityParams[0]*exp(-fNonLinearityParams[1]/energy))+
675 ((fNonLinearityParams[2]/(fNonLinearityParams[3]*2.*TMath::Pi())*
676 exp(-(energy-fNonLinearityParams[4])*(energy-fNonLinearityParams[4])/(2.*fNonLinearityParams[3]*fNonLinearityParams[3]))));
682 //Non-Linearity correction (from MC with function [0]/((x+[1])^[2]))+1;
683 //fNonLinearityParams[0] = 3.11111e-02;
684 //fNonLinearityParams[1] =-5.71666e-02;
685 //fNonLinearityParams[2] = 5.67995e-01;
687 energy *= fNonLinearityParams[0]/TMath::Power(energy+fNonLinearityParams[1],fNonLinearityParams[2])+1;
693 //Same as beam test corrected, change parameters
694 //fNonLinearityParams[0] = 9.81039e-01
695 //fNonLinearityParams[1] = 1.13508e-01;
696 //fNonLinearityParams[2] = 1.00173e+00;
697 //fNonLinearityParams[3] = 9.67998e-02;
698 //fNonLinearityParams[4] = 2.19381e+02;
699 //fNonLinearityParams[5] = 6.31604e+01;
700 //fNonLinearityParams[6] = 1;
701 energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
709 //Non-Linearity correction (from Olga Data with function p0+p1*exp(-p2*E))
710 //fNonLinearityParams[0] = 1.04;
711 //fNonLinearityParams[1] = -0.1445;
712 //fNonLinearityParams[2] = 1.046;
713 energy /= (fNonLinearityParams[0]+fNonLinearityParams[1]*exp(-fNonLinearityParams[2]*energy)); //Olga function
717 case kPi0GammaConversion:
719 //Non-Linearity correction (Nicolas from Dimitri Data with function C*[1-a*exp(-b*E)])
720 //fNonLinearityParams[0] = 0.139393/0.1349766;
721 //fNonLinearityParams[1] = 0.0566186;
722 //fNonLinearityParams[2] = 0.982133;
723 energy /= fNonLinearityParams[0]*(1-fNonLinearityParams[1]*exp(-fNonLinearityParams[2]*energy));
730 //From beam test, Alexei's results, for different ZS thresholds
731 // th=30 MeV; th = 45 MeV; th = 75 MeV
732 //fNonLinearityParams[0] = 1.007; 1.003; 1.002
733 //fNonLinearityParams[1] = 0.894; 0.719; 0.797
734 //fNonLinearityParams[2] = 0.246; 0.334; 0.358
735 //Rescale the param[0] with 1.03
736 energy /= fNonLinearityParams[0]/(1+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]));
741 case kBeamTestCorrected:
743 //From beam test, corrected for material between beam and EMCAL
744 //fNonLinearityParams[0] = 0.99078
745 //fNonLinearityParams[1] = 0.161499;
746 //fNonLinearityParams[2] = 0.655166;
747 //fNonLinearityParams[3] = 0.134101;
748 //fNonLinearityParams[4] = 163.282;
749 //fNonLinearityParams[5] = 23.6904;
750 //fNonLinearityParams[6] = 0.978;
751 energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
757 AliDebug(2,"No correction on the energy\n");
765 //__________________________________________________
766 void AliEMCALRecoUtils::InitNonLinearityParam()
768 //Initialising Non Linearity Parameters
770 if(fNonLinearityFunction == kPi0MC)
772 fNonLinearityParams[0] = 1.014;
773 fNonLinearityParams[1] = -0.03329;
774 fNonLinearityParams[2] = -0.3853;
775 fNonLinearityParams[3] = 0.5423;
776 fNonLinearityParams[4] = -0.4335;
779 if(fNonLinearityFunction == kPi0MCv2)
781 fNonLinearityParams[0] = 3.11111e-02;
782 fNonLinearityParams[1] =-5.71666e-02;
783 fNonLinearityParams[2] = 5.67995e-01;
786 if(fNonLinearityFunction == kPi0MCv3)
788 fNonLinearityParams[0] = 9.81039e-01;
789 fNonLinearityParams[1] = 1.13508e-01;
790 fNonLinearityParams[2] = 1.00173e+00;
791 fNonLinearityParams[3] = 9.67998e-02;
792 fNonLinearityParams[4] = 2.19381e+02;
793 fNonLinearityParams[5] = 6.31604e+01;
794 fNonLinearityParams[6] = 1;
797 if(fNonLinearityFunction == kPi0GammaGamma)
799 fNonLinearityParams[0] = 1.04;
800 fNonLinearityParams[1] = -0.1445;
801 fNonLinearityParams[2] = 1.046;
804 if(fNonLinearityFunction == kPi0GammaConversion)
806 fNonLinearityParams[0] = 0.139393;
807 fNonLinearityParams[1] = 0.0566186;
808 fNonLinearityParams[2] = 0.982133;
811 if(fNonLinearityFunction == kBeamTest)
813 if(fNonLinearThreshold == 30)
815 fNonLinearityParams[0] = 1.007;
816 fNonLinearityParams[1] = 0.894;
817 fNonLinearityParams[2] = 0.246;
819 if(fNonLinearThreshold == 45)
821 fNonLinearityParams[0] = 1.003;
822 fNonLinearityParams[1] = 0.719;
823 fNonLinearityParams[2] = 0.334;
825 if(fNonLinearThreshold == 75)
827 fNonLinearityParams[0] = 1.002;
828 fNonLinearityParams[1] = 0.797;
829 fNonLinearityParams[2] = 0.358;
833 if(fNonLinearityFunction == kBeamTestCorrected)
835 fNonLinearityParams[0] = 0.99078;
836 fNonLinearityParams[1] = 0.161499;
837 fNonLinearityParams[2] = 0.655166;
838 fNonLinearityParams[3] = 0.134101;
839 fNonLinearityParams[4] = 163.282;
840 fNonLinearityParams[5] = 23.6904;
841 fNonLinearityParams[6] = 0.978;
845 //_________________________________________________________
846 Float_t AliEMCALRecoUtils::GetDepth(const Float_t energy,
847 const Int_t iParticle,
848 const Int_t iSM) const
850 //Calculate shower depth for a given cluster energy and particle type
856 Float_t arg = energy*1000/ ecr; //Multiply energy by 1000 to transform to MeV
864 depth = x0 * (TMath::Log(arg) + 0.5);
871 depth = x0 * (TMath::Log(arg) - 0.5);
879 gGeoManager->cd("ALIC_1/XEN1_1");
880 TGeoNode *geoXEn1 = gGeoManager->GetCurrentNode();
881 TGeoNodeMatrix *geoSM = dynamic_cast<TGeoNodeMatrix *>(geoXEn1->GetDaughter(iSM));
884 TGeoVolume *geoSMVol = geoSM->GetVolume();
885 TGeoShape *geoSMShape = geoSMVol->GetShape();
886 TGeoBBox *geoBox = dynamic_cast<TGeoBBox *>(geoSMShape);
887 if(geoBox) depth = 0.5 * geoBox->GetDX()*2 ;
888 else AliFatal("Null GEANT box");
890 else AliFatal("NULL GEANT node matrix");
897 depth = x0 * (TMath::Log(arg) - 0.5);
906 depth = x0 * (TMath::Log(arg) + 0.5);
912 //____________________________________________________________________
913 void AliEMCALRecoUtils::GetMaxEnergyCell(const AliEMCALGeometry *geom,
914 AliVCaloCells* cells,
915 const AliVCluster* clu,
922 //For a given CaloCluster gets the absId of the cell
923 //with maximum energy deposit.
926 Double_t eCell = -1.;
927 Float_t fraction = 1.;
928 Float_t recalFactor = 1.;
929 Int_t cellAbsId = -1 ;
938 AliInfo("Cluster pointer null!");
939 absId=-1; iSupMod0=-1, ieta = -1; iphi = -1; shared = -1;
943 for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++)
945 cellAbsId = clu->GetCellAbsId(iDig);
946 fraction = clu->GetCellAmplitudeFraction(iDig);
947 //printf("a Cell %d, id, %d, amp %f, fraction %f\n",iDig,cellAbsId,cells->GetCellAmplitude(cellAbsId),fraction);
948 if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
949 geom->GetCellIndex(cellAbsId,iSupMod,iTower,iIphi,iIeta);
950 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);
955 else if(iSupMod0!=iSupMod)
958 //printf("AliEMCALRecoUtils::GetMaxEnergyCell() - SHARED CLUSTER\n");
960 if(!fCellsRecalibrated && IsRecalibrationOn())
962 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
964 eCell = cells->GetCellAmplitude(cellAbsId)*fraction*recalFactor;
965 //printf("b Cell %d, id, %d, amp %f, fraction %f\n",iDig,cellAbsId,eCell,fraction);
970 //printf("\t new max: cell %d, e %f, ecell %f\n",maxId, eMax,eCell);
974 //Get from the absid the supermodule, tower and eta/phi numbers
975 geom->GetCellIndex(absId,iSupMod,iTower,iIphi,iIeta);
976 //Gives SuperModule and Tower numbers
977 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,
978 iIphi, iIeta,iphi,ieta);
979 //printf("Max id %d, iSM %d, col %d, row %d\n",absId,iSupMod,ieta,iphi);
980 //printf("Max end---\n");
983 //______________________________________
984 void AliEMCALRecoUtils::InitParameters()
986 // Initialize data members with default values
988 fParticleType = kPhoton;
989 fPosAlgo = kUnchanged;
992 fNonLinearityFunction = kNoCorrection;
993 fNonLinearThreshold = 30;
995 fExoticCellFraction = 0.97;
996 fExoticCellDiffTime = 1e6;
997 fExoticCellMinAmplitude = 0.5;
1001 fCutEtaPhiSum = kTRUE;
1002 fCutEtaPhiSeparate = kFALSE;
1008 fClusterWindow = 100;
1013 fTrackCutsType = kLooseCut;
1016 fCutMinNClusterTPC = -1;
1017 fCutMinNClusterITS = -1;
1019 fCutMaxChi2PerClusterTPC = 1e10;
1020 fCutMaxChi2PerClusterITS = 1e10;
1022 fCutRequireTPCRefit = kFALSE;
1023 fCutRequireITSRefit = kFALSE;
1024 fCutAcceptKinkDaughters = kFALSE;
1026 fCutMaxDCAToVertexXY = 1e10;
1027 fCutMaxDCAToVertexZ = 1e10;
1028 fCutDCAToVertex2D = kFALSE;
1030 fCutRequireITSStandAlone = kFALSE; //MARCEL
1031 fCutRequireITSpureSA = kFALSE; //Marcel
1033 //Misalignment matrices
1034 for(Int_t i = 0; i < 15 ; i++)
1036 fMisalTransShift[i] = 0.;
1037 fMisalRotShift[i] = 0.;
1041 for(Int_t i = 0; i < 7 ; i++) fNonLinearityParams[i] = 0.;
1043 //For kBeamTestCorrected case, but default is no correction
1044 fNonLinearityParams[0] = 0.99078;
1045 fNonLinearityParams[1] = 0.161499;
1046 fNonLinearityParams[2] = 0.655166;
1047 fNonLinearityParams[3] = 0.134101;
1048 fNonLinearityParams[4] = 163.282;
1049 fNonLinearityParams[5] = 23.6904;
1050 fNonLinearityParams[6] = 0.978;
1052 //For kPi0GammaGamma case
1053 //fNonLinearityParams[0] = 0.1457/0.1349766/1.038;
1054 //fNonLinearityParams[1] = -0.02024/0.1349766/1.038;
1055 //fNonLinearityParams[2] = 1.046;
1057 //Cluster energy smearing
1058 fSmearClusterEnergy = kFALSE;
1059 fSmearClusterParam[0] = 0.07; // * sqrt E term
1060 fSmearClusterParam[1] = 0.00; // * E term
1061 fSmearClusterParam[2] = 0.00; // constant
1064 //_____________________________________________________
1065 void AliEMCALRecoUtils::InitEMCALRecalibrationFactors()
1067 //Init EMCAL recalibration factors
1068 AliDebug(2,"AliCalorimeterUtils::InitEMCALRecalibrationFactors()");
1069 //In order to avoid rewriting the same histograms
1070 Bool_t oldStatus = TH1::AddDirectoryStatus();
1071 TH1::AddDirectory(kFALSE);
1073 fEMCALRecalibrationFactors = new TObjArray(12);
1074 for (int i = 0; i < 12; i++)
1075 fEMCALRecalibrationFactors->Add(new TH2F(Form("EMCALRecalFactors_SM%d",i),
1076 Form("EMCALRecalFactors_SM%d",i), 48, 0, 48, 24, 0, 24));
1077 //Init the histograms with 1
1078 for (Int_t sm = 0; sm < 12; sm++)
1080 for (Int_t i = 0; i < 48; i++)
1082 for (Int_t j = 0; j < 24; j++)
1084 SetEMCALChannelRecalibrationFactor(sm,i,j,1.);
1089 fEMCALRecalibrationFactors->SetOwner(kTRUE);
1090 fEMCALRecalibrationFactors->Compress();
1092 //In order to avoid rewriting the same histograms
1093 TH1::AddDirectory(oldStatus);
1096 //_________________________________________________________
1097 void AliEMCALRecoUtils::InitEMCALTimeRecalibrationFactors()
1099 //Init EMCAL recalibration factors
1100 AliDebug(2,"AliCalorimeterUtils::InitEMCALRecalibrationFactors()");
1101 //In order to avoid rewriting the same histograms
1102 Bool_t oldStatus = TH1::AddDirectoryStatus();
1103 TH1::AddDirectory(kFALSE);
1105 fEMCALTimeRecalibrationFactors = new TObjArray(4);
1106 for (int i = 0; i < 4; i++)
1107 fEMCALTimeRecalibrationFactors->Add(new TH1F(Form("hAllTimeAvBC%d",i),
1108 Form("hAllTimeAvBC%d",i),
1109 48*24*12,0.,48*24*12) );
1110 //Init the histograms with 1
1111 for (Int_t bc = 0; bc < 4; bc++)
1113 for (Int_t i = 0; i < 48*24*12; i++)
1114 SetEMCALChannelTimeRecalibrationFactor(bc,i,0.);
1117 fEMCALTimeRecalibrationFactors->SetOwner(kTRUE);
1118 fEMCALTimeRecalibrationFactors->Compress();
1120 //In order to avoid rewriting the same histograms
1121 TH1::AddDirectory(oldStatus);
1124 //____________________________________________________
1125 void AliEMCALRecoUtils::InitEMCALBadChannelStatusMap()
1127 //Init EMCAL bad channels map
1128 AliDebug(2,"AliEMCALRecoUtils::InitEMCALBadChannelStatusMap()");
1129 //In order to avoid rewriting the same histograms
1130 Bool_t oldStatus = TH1::AddDirectoryStatus();
1131 TH1::AddDirectory(kFALSE);
1133 fEMCALBadChannelMap = new TObjArray(12);
1134 //TH2F * hTemp = new TH2I("EMCALBadChannelMap","EMCAL SuperModule bad channel map", 48, 0, 48, 24, 0, 24);
1135 for (int i = 0; i < 12; i++)
1137 fEMCALBadChannelMap->Add(new TH2I(Form("EMCALBadChannelMap_Mod%d",i),Form("EMCALBadChannelMap_Mod%d",i), 48, 0, 48, 24, 0, 24));
1140 fEMCALBadChannelMap->SetOwner(kTRUE);
1141 fEMCALBadChannelMap->Compress();
1143 //In order to avoid rewriting the same histograms
1144 TH1::AddDirectory(oldStatus);
1147 //____________________________________________________________________________
1148 void AliEMCALRecoUtils::RecalibrateClusterEnergy(const AliEMCALGeometry* geom,
1149 AliVCluster * cluster,
1150 AliVCaloCells * cells,
1153 // Recalibrate the cluster energy and Time, considering the recalibration map
1154 // and the energy of the cells and time that compose the cluster.
1155 // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1159 AliInfo("Cluster pointer null!");
1163 //Get the cluster number of cells and list of absId, check what kind of cluster do we have.
1164 UShort_t * index = cluster->GetCellsAbsId() ;
1165 Double_t * fraction = cluster->GetCellsAmplitudeFraction() ;
1166 Int_t ncells = cluster->GetNCells();
1168 //Initialize some used variables
1171 Int_t icol =-1, irow =-1, imod=1;
1172 Float_t factor = 1, frac = 0;
1173 Int_t absIdMax = -1;
1176 //Loop on the cells, get the cell amplitude and recalibration factor, multiply and and to the new energy
1177 for(Int_t icell = 0; icell < ncells; icell++)
1179 absId = index[icell];
1180 frac = fraction[icell];
1181 if(frac < 1e-5) frac = 1; //in case of EMCAL, this is set as 0 since unfolding is off
1183 if(!fCellsRecalibrated && IsRecalibrationOn())
1186 Int_t iTower = -1, iIphi = -1, iIeta = -1;
1187 geom->GetCellIndex(absId,imod,iTower,iIphi,iIeta);
1188 if(fEMCALRecalibrationFactors->GetEntries() <= imod) continue;
1189 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,irow,icol);
1190 factor = GetEMCALChannelRecalibrationFactor(imod,icol,irow);
1192 AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - recalibrate cell: module %d, col %d, row %d, cell fraction %f,recalibration factor %f, cell energy %f\n",
1193 imod,icol,irow,frac,factor,cells->GetCellAmplitude(absId)));
1197 energy += cells->GetCellAmplitude(absId)*factor*frac;
1199 if(emax < cells->GetCellAmplitude(absId)*factor*frac)
1201 emax = cells->GetCellAmplitude(absId)*factor*frac;
1206 AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - Energy before %f, after %f \n",cluster->E(),energy));
1208 cluster->SetE(energy);
1210 // Recalculate time of cluster
1211 Double_t timeorg = cluster->GetTOF();
1213 Double_t time = cells->GetCellTime(absIdMax);
1214 if(!fCellsRecalibrated && IsTimeRecalibrationOn())
1215 RecalibrateCellTime(absIdMax,bc,time);
1217 cluster->SetTOF(time);
1219 AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - Time before %f, after %f \n",timeorg,cluster->GetTOF()));
1222 //_____________________________________________________________
1223 void AliEMCALRecoUtils::RecalibrateCells(AliVCaloCells * cells,
1226 // Recalibrate the cells time and energy, considering the recalibration map and the energy
1227 // of the cells that compose the cluster.
1228 // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1230 if(!IsRecalibrationOn() && !IsTimeRecalibrationOn() && !IsBadChannelsRemovalSwitchedOn()) return;
1234 AliInfo("Cells pointer null!");
1239 Bool_t accept = kFALSE;
1242 Double_t ecellin = 0;
1243 Double_t tcellin = 0;
1244 Short_t mclabel = -1;
1247 Int_t nEMcell = cells->GetNumberOfCells() ;
1248 for (Int_t iCell = 0; iCell < nEMcell; iCell++)
1250 cells->GetCell( iCell, absId, ecellin, tcellin, mclabel, efrac );
1252 accept = AcceptCalibrateCell(absId, bc, ecell ,tcell ,cells);
1260 cells->SetCell(iCell,absId,ecell, tcell, mclabel, efrac);
1263 fCellsRecalibrated = kTRUE;
1266 //_______________________________________________________________________________________________________
1267 void AliEMCALRecoUtils::RecalibrateCellTime(const Int_t absId, const Int_t bc, Double_t & celltime) const
1269 // Recalibrate time of cell with absID considering the recalibration map
1270 // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1272 if(!fCellsRecalibrated && IsTimeRecalibrationOn() && bc >= 0)
1274 celltime -= GetEMCALChannelTimeRecalibrationFactor(bc%4,absId)*1.e-9; ;
1278 //______________________________________________________________________________
1279 void AliEMCALRecoUtils::RecalculateClusterPosition(const AliEMCALGeometry *geom,
1280 AliVCaloCells* cells,
1283 //For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1287 AliInfo("Cluster pointer null!");
1291 if (fPosAlgo==kPosTowerGlobal) RecalculateClusterPositionFromTowerGlobal( geom, cells, clu);
1292 else if(fPosAlgo==kPosTowerIndex) RecalculateClusterPositionFromTowerIndex ( geom, cells, clu);
1293 else AliDebug(2,"Algorithm to recalculate position not selected, do nothing.");
1296 //_____________________________________________________________________________________________
1297 void AliEMCALRecoUtils::RecalculateClusterPositionFromTowerGlobal(const AliEMCALGeometry *geom,
1298 AliVCaloCells* cells,
1301 // For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1302 // The algorithm is a copy of what is done in AliEMCALRecPoint
1304 Double_t eCell = 0.;
1305 Float_t fraction = 1.;
1306 Float_t recalFactor = 1.;
1309 Int_t iTower = -1, iIphi = -1, iIeta = -1;
1310 Int_t iSupModMax = -1, iSM=-1, iphi = -1, ieta = -1;
1311 Float_t weight = 0., totalWeight=0.;
1312 Float_t newPos[3] = {0,0,0};
1313 Double_t pLocal[3], pGlobal[3];
1314 Bool_t shared = kFALSE;
1316 Float_t clEnergy = clu->E(); //Energy already recalibrated previously
1319 GetMaxEnergyCell(geom, cells, clu, absId, iSupModMax, ieta, iphi,shared);
1320 Double_t depth = GetDepth(clEnergy,fParticleType,iSupModMax) ;
1322 //printf("** Cluster energy %f, ncells %d, depth %f\n",clEnergy,clu->GetNCells(),depth);
1324 for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++)
1326 absId = clu->GetCellAbsId(iDig);
1327 fraction = clu->GetCellAmplitudeFraction(iDig);
1328 if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1330 if (!fCellsRecalibrated)
1332 geom->GetCellIndex(absId,iSM,iTower,iIphi,iIeta);
1333 geom->GetCellPhiEtaIndexInSModule(iSM,iTower,iIphi, iIeta,iphi,ieta);
1335 if(IsRecalibrationOn())
1337 recalFactor = GetEMCALChannelRecalibrationFactor(iSM,ieta,iphi);
1341 eCell = cells->GetCellAmplitude(absId)*fraction*recalFactor;
1343 weight = GetCellWeight(eCell,clEnergy);
1344 totalWeight += weight;
1346 geom->RelPosCellInSModule(absId,depth,pLocal[0],pLocal[1],pLocal[2]);
1347 //printf("pLocal (%f,%f,%f), SM %d, absId %d\n",pLocal[0],pLocal[1],pLocal[2],iSupModMax,absId);
1348 geom->GetGlobal(pLocal,pGlobal,iSupModMax);
1349 //printf("pLocal (%f,%f,%f)\n",pGlobal[0],pGlobal[1],pGlobal[2]);
1351 for(int i=0; i<3; i++ ) newPos[i] += (weight*pGlobal[i]);
1356 for(int i=0; i<3; i++ ) newPos[i] /= totalWeight;
1359 //Float_t pos[]={0,0,0};
1360 //clu->GetPosition(pos);
1361 //printf("OldPos : %2.3f,%2.3f,%2.3f\n",pos[0],pos[1],pos[2]);
1362 //printf("NewPos : %2.3f,%2.3f,%2.3f\n",newPos[0],newPos[1],newPos[2]);
1364 if(iSupModMax > 1) //sector 1
1366 newPos[0] +=fMisalTransShift[3];//-=3.093;
1367 newPos[1] +=fMisalTransShift[4];//+=6.82;
1368 newPos[2] +=fMisalTransShift[5];//+=1.635;
1369 //printf(" + : %2.3f,%2.3f,%2.3f\n",fMisalTransShift[3],fMisalTransShift[4],fMisalTransShift[5]);
1372 newPos[0] +=fMisalTransShift[0];//+=1.134;
1373 newPos[1] +=fMisalTransShift[1];//+=8.2;
1374 newPos[2] +=fMisalTransShift[2];//+=1.197;
1375 //printf(" + : %2.3f,%2.3f,%2.3f\n",fMisalTransShift[0],fMisalTransShift[1],fMisalTransShift[2]);
1377 //printf("NewPos : %2.3f,%2.3f,%2.3f\n",newPos[0],newPos[1],newPos[2]);
1379 clu->SetPosition(newPos);
1382 //____________________________________________________________________________________________
1383 void AliEMCALRecoUtils::RecalculateClusterPositionFromTowerIndex(const AliEMCALGeometry *geom,
1384 AliVCaloCells* cells,
1387 // For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1388 // The algorithm works with the tower indeces, averages the indeces and from them it calculates the global position
1390 Double_t eCell = 1.;
1391 Float_t fraction = 1.;
1392 Float_t recalFactor = 1.;
1396 Int_t iIphi = -1, iIeta = -1;
1397 Int_t iSupMod = -1, iSupModMax = -1;
1398 Int_t iphi = -1, ieta =-1;
1399 Bool_t shared = kFALSE;
1401 Float_t clEnergy = clu->E(); //Energy already recalibrated previously.
1404 GetMaxEnergyCell(geom, cells, clu, absId, iSupModMax, ieta, iphi,shared);
1405 Float_t depth = GetDepth(clEnergy,fParticleType,iSupMod) ;
1407 Float_t weight = 0., weightedCol = 0., weightedRow = 0., totalWeight=0.;
1408 Bool_t areInSameSM = kTRUE; //exclude clusters with cells in different SMs for now
1409 Int_t startingSM = -1;
1411 for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++)
1413 absId = clu->GetCellAbsId(iDig);
1414 fraction = clu->GetCellAmplitudeFraction(iDig);
1415 if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1417 if (iDig==0) startingSM = iSupMod;
1418 else if(iSupMod != startingSM) areInSameSM = kFALSE;
1420 eCell = cells->GetCellAmplitude(absId);
1422 geom->GetCellIndex(absId,iSupMod,iTower,iIphi,iIeta);
1423 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);
1425 if (!fCellsRecalibrated)
1427 if(IsRecalibrationOn())
1429 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1433 eCell = cells->GetCellAmplitude(absId)*fraction*recalFactor;
1435 weight = GetCellWeight(eCell,clEnergy);
1436 if(weight < 0) weight = 0;
1437 totalWeight += weight;
1438 weightedCol += ieta*weight;
1439 weightedRow += iphi*weight;
1441 //printf("Max cell? cell %d, amplitude org %f, fraction %f, recalibration %f, amplitude new %f \n",cellAbsId, cells->GetCellAmplitude(cellAbsId), fraction, recalFactor, eCell) ;
1444 Float_t xyzNew[]={0.,0.,0.};
1445 if(areInSameSM == kTRUE)
1447 //printf("In Same SM\n");
1448 weightedCol = weightedCol/totalWeight;
1449 weightedRow = weightedRow/totalWeight;
1450 geom->RecalculateTowerPosition(weightedRow, weightedCol, iSupModMax, depth, fMisalTransShift, fMisalRotShift, xyzNew);
1454 //printf("In Different SM\n");
1455 geom->RecalculateTowerPosition(iphi, ieta, iSupModMax, depth, fMisalTransShift, fMisalRotShift, xyzNew);
1458 clu->SetPosition(xyzNew);
1461 //___________________________________________________________________________________________
1462 void AliEMCALRecoUtils::RecalculateClusterDistanceToBadChannel(const AliEMCALGeometry * geom,
1463 AliVCaloCells* cells,
1464 AliVCluster * cluster)
1466 //re-evaluate distance to bad channel with updated bad map
1468 if(!fRecalDistToBadChannels) return;
1472 AliInfo("Cluster pointer null!");
1476 //Get channels map of the supermodule where the cluster is.
1477 Int_t absIdMax = -1, iSupMod =-1, icolM = -1, irowM = -1;
1478 Bool_t shared = kFALSE;
1479 GetMaxEnergyCell(geom, cells, cluster, absIdMax, iSupMod, icolM, irowM, shared);
1480 TH2D* hMap = (TH2D*)fEMCALBadChannelMap->At(iSupMod);
1483 Float_t minDist = 10000.;
1486 //Loop on tower status map
1487 for(Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++)
1489 for(Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++)
1491 //Check if tower is bad.
1492 if(hMap->GetBinContent(icol,irow)==0) continue;
1493 //printf("AliEMCALRecoUtils::RecalculateDistanceToBadChannels() - \n \t Bad channel in SM %d, col %d, row %d, \n \t Cluster max in col %d, row %d\n",
1494 // iSupMod,icol, irow, icolM,irowM);
1496 dRrow=TMath::Abs(irowM-irow);
1497 dRcol=TMath::Abs(icolM-icol);
1498 dist=TMath::Sqrt(dRrow*dRrow+dRcol*dRcol);
1501 //printf("MIN DISTANCE TO BAD %2.2f\n",dist);
1507 //In case the cluster is shared by 2 SuperModules, need to check the map of the second Super Module
1511 Int_t iSupMod2 = -1;
1513 //The only possible combinations are (0,1), (2,3) ... (8,9)
1514 if(iSupMod%2) iSupMod2 = iSupMod-1;
1515 else iSupMod2 = iSupMod+1;
1516 hMap2 = (TH2D*)fEMCALBadChannelMap->At(iSupMod2);
1518 //Loop on tower status map of second super module
1519 for(Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++)
1521 for(Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++)
1523 //Check if tower is bad.
1524 if(hMap2->GetBinContent(icol,irow)==0) continue;
1525 //printf("AliEMCALRecoUtils::RecalculateDistanceToBadChannels(shared) - \n \t Bad channel in SM %d, col %d, row %d \n \t Cluster max in SM %d, col %d, row %d\n",
1526 // iSupMod2,icol, irow,iSupMod,icolM,irowM);
1527 dRrow=TMath::Abs(irow-irowM);
1531 dRcol=TMath::Abs(icol-(AliEMCALGeoParams::fgkEMCALCols+icolM));
1534 dRcol=TMath::Abs(AliEMCALGeoParams::fgkEMCALCols+icol-icolM);
1537 dist=TMath::Sqrt(dRrow*dRrow+dRcol*dRcol);
1538 if(dist < minDist) minDist = dist;
1541 }// shared cluster in 2 SuperModules
1543 AliDebug(2,Form("Max cluster cell (SM,col,row)=(%d %d %d) - Distance to Bad Channel %2.2f",iSupMod, icolM, irowM, minDist));
1544 cluster->SetDistanceToBadChannel(minDist);
1547 //__________________________________________________________________
1548 void AliEMCALRecoUtils::RecalculateClusterPID(AliVCluster * cluster)
1550 //re-evaluate identification parameters with bayesian
1554 AliInfo("Cluster pointer null!");
1558 if ( cluster->GetM02() != 0)
1559 fPIDUtils->ComputePID(cluster->E(),cluster->GetM02());
1561 Float_t pidlist[AliPID::kSPECIESCN+1];
1562 for(Int_t i = 0; i < AliPID::kSPECIESCN+1; i++) pidlist[i] = fPIDUtils->GetPIDFinal(i);
1564 cluster->SetPID(pidlist);
1567 //___________________________________________________________________________________________________________________
1568 void AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(const AliEMCALGeometry * geom,
1569 AliVCaloCells* cells,
1570 AliVCluster * cluster,
1571 Float_t & l0, Float_t & l1,
1572 Float_t & disp, Float_t & dEta, Float_t & dPhi,
1573 Float_t & sEta, Float_t & sPhi, Float_t & sEtaPhi)
1575 // Calculates new center of gravity in the local EMCAL-module coordinates
1576 // and tranfers into global ALICE coordinates
1577 // Calculates Dispersion and main axis
1581 AliInfo("Cluster pointer null!");
1585 Double_t eCell = 0.;
1586 Float_t fraction = 1.;
1587 Float_t recalFactor = 1.;
1595 Double_t etai = -1.;
1596 Double_t phii = -1.;
1601 Double_t etaMean = 0.;
1602 Double_t phiMean = 0.;
1605 for(Int_t iDigit=0; iDigit < cluster->GetNCells(); iDigit++)
1607 //Get from the absid the supermodule, tower and eta/phi numbers
1608 geom->GetCellIndex(cluster->GetCellAbsId(iDigit),iSupMod,iTower,iIphi,iIeta);
1609 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi,iIeta, iphi,ieta);
1611 //Get the cell energy, if recalibration is on, apply factors
1612 fraction = cluster->GetCellAmplitudeFraction(iDigit);
1613 if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1615 if (!fCellsRecalibrated)
1617 if(IsRecalibrationOn())
1619 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1623 eCell = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1625 if(cluster->E() > 0 && eCell > 0)
1627 w = GetCellWeight(eCell,cluster->E());
1629 etai=(Double_t)ieta;
1630 phii=(Double_t)iphi;
1637 sEta += w * etai * etai ;
1638 etaMean += w * etai ;
1639 sPhi += w * phii * phii ;
1640 phiMean += w * phii ;
1641 sEtaPhi += w * etai * phii ;
1645 AliError(Form("Wrong energy %f and/or amplitude %f\n", eCell, cluster->E()));
1648 //Normalize to the weight
1655 AliError(Form("Wrong weight %f\n", wtot));
1657 //Calculate dispersion
1658 for(Int_t iDigit=0; iDigit < cluster->GetNCells(); iDigit++)
1660 //Get from the absid the supermodule, tower and eta/phi numbers
1661 geom->GetCellIndex(cluster->GetCellAbsId(iDigit),iSupMod,iTower,iIphi,iIeta);
1662 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi,iIeta, iphi,ieta);
1664 //Get the cell energy, if recalibration is on, apply factors
1665 fraction = cluster->GetCellAmplitudeFraction(iDigit);
1666 if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1667 if (IsRecalibrationOn())
1669 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1671 eCell = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1673 if(cluster->E() > 0 && eCell > 0)
1675 w = GetCellWeight(eCell,cluster->E());
1677 etai=(Double_t)ieta;
1678 phii=(Double_t)iphi;
1681 disp += w *((etai-etaMean)*(etai-etaMean)+(phii-phiMean)*(phii-phiMean));
1682 dEta += w * (etai-etaMean)*(etai-etaMean) ;
1683 dPhi += w * (phii-phiMean)*(phii-phiMean) ;
1687 AliError(Form("Wrong energy %f and/or amplitude %f\n", eCell, cluster->E()));
1690 //Normalize to the weigth and set shower shape parameters
1691 if (wtot > 0 && nstat > 1)
1700 sEta -= etaMean * etaMean ;
1701 sPhi -= phiMean * phiMean ;
1702 sEtaPhi -= etaMean * phiMean ;
1704 l0 = (0.5 * (sEta + sPhi) + TMath::Sqrt( 0.25 * (sEta - sPhi) * (sEta - sPhi) + sEtaPhi * sEtaPhi ));
1705 l1 = (0.5 * (sEta + sPhi) - TMath::Sqrt( 0.25 * (sEta - sPhi) * (sEta - sPhi) + sEtaPhi * sEtaPhi ));
1711 dEta = 0. ; dPhi = 0. ; disp = 0. ;
1712 sEta = 0. ; sPhi = 0. ; sEtaPhi = 0. ;
1717 //____________________________________________________________________________________________
1718 void AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(const AliEMCALGeometry * geom,
1719 AliVCaloCells* cells,
1720 AliVCluster * cluster)
1722 // Calculates new center of gravity in the local EMCAL-module coordinates
1723 // and tranfers into global ALICE coordinates
1724 // Calculates Dispersion and main axis and puts them into the cluster
1726 Float_t l0 = 0., l1 = 0.;
1727 Float_t disp = 0., dEta = 0., dPhi = 0.;
1728 Float_t sEta = 0., sPhi = 0., sEtaPhi = 0.;
1730 AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(geom,cells,cluster,l0,l1,disp,
1731 dEta, dPhi, sEta, sPhi, sEtaPhi);
1733 cluster->SetM02(l0);
1734 cluster->SetM20(l1);
1735 if(disp > 0. ) cluster->SetDispersion(TMath::Sqrt(disp)) ;
1739 //____________________________________________________________________________
1740 void AliEMCALRecoUtils::FindMatches(AliVEvent *event,
1741 TObjArray * clusterArr,
1742 const AliEMCALGeometry *geom)
1744 //This function should be called before the cluster loop
1745 //Before call this function, please recalculate the cluster positions
1746 //Given the input event, loop over all the tracks, select the closest cluster as matched with fCutR
1747 //Store matched cluster indexes and residuals
1749 fMatchedTrackIndex ->Reset();
1750 fMatchedClusterIndex->Reset();
1751 fResidualPhi->Reset();
1752 fResidualEta->Reset();
1754 fMatchedTrackIndex ->Set(1000);
1755 fMatchedClusterIndex->Set(1000);
1756 fResidualPhi->Set(1000);
1757 fResidualEta->Set(1000);
1759 AliESDEvent* esdevent = dynamic_cast<AliESDEvent*> (event);
1760 AliAODEvent* aodevent = dynamic_cast<AliAODEvent*> (event);
1762 // init the magnetic field if not already on
1763 if(!TGeoGlobalMagField::Instance()->GetField())
1765 AliInfo("Init the magnetic field\n");
1768 esdevent->InitMagneticField();
1772 Double_t curSol = 30000*aodevent->GetMagneticField()/5.00668;
1773 Double_t curDip = 6000 *aodevent->GetMuonMagFieldScale();
1774 AliMagF *field = AliMagF::CreateFieldMap(curSol,curDip);
1775 TGeoGlobalMagField::Instance()->SetField(field);
1779 AliInfo("Mag Field not initialized, null esd/aod evetn pointers");
1785 UInt_t mask1 = esdevent->GetESDRun()->GetDetectorsInDAQ();
1786 UInt_t mask2 = esdevent->GetESDRun()->GetDetectorsInReco();
1787 Bool_t desc1 = (mask1 >> 3) & 0x1;
1788 Bool_t desc2 = (mask2 >> 3) & 0x1;
1789 if (desc1==0 || desc2==0) {
1790 AliError(Form("TPC not in DAQ/RECO: %u (%u)/%u (%u)",
1791 mask1, esdevent->GetESDRun()->GetDetectorsInReco(),
1792 mask2, esdevent->GetESDRun()->GetDetectorsInDAQ()));
1797 TObjArray *clusterArray = 0x0;
1800 clusterArray = new TObjArray(event->GetNumberOfCaloClusters());
1801 for(Int_t icl=0; icl<event->GetNumberOfCaloClusters(); icl++)
1803 AliVCluster *cluster = (AliVCluster*) event->GetCaloCluster(icl);
1804 if(geom && !IsGoodCluster(cluster,geom,(AliVCaloCells*)event->GetEMCALCells())) continue;
1805 clusterArray->AddAt(cluster,icl);
1811 for (Int_t i=0; i<21;i++) cv[i]=0;
1812 for(Int_t itr=0; itr<event->GetNumberOfTracks(); itr++)
1814 AliExternalTrackParam *trackParam = 0;
1816 //If the input event is ESD, the starting point for extrapolation is TPCOut, if available, or TPCInner
1817 AliESDtrack *esdTrack = 0;
1818 AliAODTrack *aodTrack = 0;
1821 esdTrack = esdevent->GetTrack(itr);
1822 if(!esdTrack) continue;
1823 if(!IsAccepted(esdTrack)) continue;
1824 if(esdTrack->Pt()<fCutMinTrackPt) continue;
1825 Double_t phi = esdTrack->Phi()*TMath::RadToDeg();
1826 if(TMath::Abs(esdTrack->Eta())>0.8 || phi <= 20 || phi >= 240 ) continue;
1828 trackParam = const_cast<AliExternalTrackParam*>(esdTrack->GetInnerParam()); // if TPC Available
1830 trackParam = new AliExternalTrackParam(*esdTrack); // If ITS Track Standing alone
1833 //If the input event is AOD, the starting point for extrapolation is at vertex
1834 //AOD tracks are selected according to its filterbit.
1837 aodTrack = aodevent->GetTrack(itr);
1838 if(!aodTrack) continue;
1839 if(!aodTrack->TestFilterMask(fAODFilterMask)) continue; //Select AOD tracks that fulfill GetStandardITSTPCTrackCuts2010()
1840 if(aodTrack->Pt()<fCutMinTrackPt) continue;
1841 Double_t phi = aodTrack->Phi()*TMath::RadToDeg();
1842 if(TMath::Abs(aodTrack->Eta())>0.8 || phi <= 20 || phi >= 240 ) continue;
1843 Double_t pos[3],mom[3];
1844 aodTrack->GetXYZ(pos);
1845 aodTrack->GetPxPyPz(mom);
1846 AliDebug(5,Form("aod track: i=%d | pos=(%5.4f,%5.4f,%5.4f) | mom=(%5.4f,%5.4f,%5.4f) | charge=%d\n",itr,pos[0],pos[1],pos[2],mom[0],mom[1],mom[2],aodTrack->Charge()));
1847 trackParam= new AliExternalTrackParam(pos,mom,cv,aodTrack->Charge());
1850 //Return if the input data is not "AOD" or "ESD"
1853 printf("Wrong input data type! Should be \"AOD\" or \"ESD\"\n");
1856 clusterArray->Clear();
1857 delete clusterArray;
1862 if(!trackParam) continue;
1864 //Extrapolate the track to EMCal surface
1865 AliExternalTrackParam emcalParam(*trackParam);
1867 if(!ExtrapolateTrackToEMCalSurface(&emcalParam, 430., fMass, fStepSurface, eta, phi))
1869 if(aodevent && trackParam) delete trackParam;
1875 // esdTrack->SetOuterParam(&emcalParam,AliExternalTrackParam::kMultSec);
1878 if(TMath::Abs(eta)>0.75 || (phi) < 70*TMath::DegToRad() || (phi) > 190*TMath::DegToRad())
1880 if(aodevent && trackParam) delete trackParam;
1885 //Find matched clusters
1887 Float_t dEta = -999, dPhi = -999;
1890 index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArray, dEta, dPhi);
1894 index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArr, dEta, dPhi);
1899 fMatchedTrackIndex ->AddAt(itr,matched);
1900 fMatchedClusterIndex ->AddAt(index,matched);
1901 fResidualEta ->AddAt(dEta,matched);
1902 fResidualPhi ->AddAt(dPhi,matched);
1905 if(aodevent && trackParam) delete trackParam;
1910 clusterArray->Clear();
1911 delete clusterArray;
1914 AliDebug(2,Form("Number of matched pairs = %d !\n",matched));
1916 fMatchedTrackIndex ->Set(matched);
1917 fMatchedClusterIndex ->Set(matched);
1918 fResidualPhi ->Set(matched);
1919 fResidualEta ->Set(matched);
1922 //________________________________________________________________________________
1923 Int_t AliEMCALRecoUtils::FindMatchedClusterInEvent(const AliESDtrack *track,
1924 const AliVEvent *event,
1925 const AliEMCALGeometry *geom,
1926 Float_t &dEta, Float_t &dPhi)
1929 // This function returns the index of matched cluster to input track
1930 // Returns -1 if no match is found
1932 Double_t phiV = track->Phi()*TMath::RadToDeg();
1933 if(TMath::Abs(track->Eta())>0.8 || phiV <= 20 || phiV >= 240 ) return index;
1934 AliExternalTrackParam *trackParam = 0;
1936 trackParam = const_cast<AliExternalTrackParam*>(track->GetInnerParam()); // If TPC
1938 trackParam = new AliExternalTrackParam(*track);
1940 if(!trackParam) return index;
1941 AliExternalTrackParam emcalParam(*trackParam);
1943 if(!ExtrapolateTrackToEMCalSurface(&emcalParam, 430., fMass, fStepSurface, eta, phi)) return index;
1944 if(TMath::Abs(eta)>0.75 || (phi) < 70*TMath::DegToRad() || (phi) > 190*TMath::DegToRad()) return index;
1946 TObjArray *clusterArr = new TObjArray(event->GetNumberOfCaloClusters());
1948 for(Int_t icl=0; icl<event->GetNumberOfCaloClusters(); icl++)
1950 AliVCluster *cluster = (AliVCluster*) event->GetCaloCluster(icl);
1951 if(geom && !IsGoodCluster(cluster,geom,(AliVCaloCells*)event->GetEMCALCells())) continue;
1952 clusterArr->AddAt(cluster,icl);
1955 index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArr, dEta, dPhi);
1956 clusterArr->Clear();
1962 //_______________________________________________________________________________________________
1963 Int_t AliEMCALRecoUtils::FindMatchedClusterInClusterArr(const AliExternalTrackParam *emcalParam,
1964 AliExternalTrackParam *trkParam,
1965 const TObjArray * clusterArr,
1966 Float_t &dEta, Float_t &dPhi)
1968 // Find matched cluster in array
1970 dEta=-999, dPhi=-999;
1971 Float_t dRMax = fCutR, dEtaMax=fCutEta, dPhiMax=fCutPhi;
1973 Float_t tmpEta=-999, tmpPhi=-999;
1975 Double_t exPos[3] = {0.,0.,0.};
1976 if(!emcalParam->GetXYZ(exPos)) return index;
1978 Float_t clsPos[3] = {0.,0.,0.};
1979 for(Int_t icl=0; icl<clusterArr->GetEntriesFast(); icl++)
1981 AliVCluster *cluster = dynamic_cast<AliVCluster*> (clusterArr->At(icl)) ;
1982 if(!cluster || !cluster->IsEMCAL()) continue;
1983 cluster->GetPosition(clsPos);
1984 Double_t dR = TMath::Sqrt(TMath::Power(exPos[0]-clsPos[0],2)+TMath::Power(exPos[1]-clsPos[1],2)+TMath::Power(exPos[2]-clsPos[2],2));
1985 if(dR > fClusterWindow) continue;
1987 AliExternalTrackParam trkPamTmp (*trkParam);//Retrieve the starting point every time before the extrapolation
1988 if(!ExtrapolateTrackToCluster(&trkPamTmp, cluster, fMass, fStepCluster, tmpEta, tmpPhi)) continue;
1991 Float_t tmpR=TMath::Sqrt(tmpEta*tmpEta + tmpPhi*tmpPhi);
2000 else if(fCutEtaPhiSeparate)
2002 if(TMath::Abs(tmpEta)<TMath::Abs(dEtaMax) && TMath::Abs(tmpPhi)<TMath::Abs(dPhiMax))
2011 printf("Error: please specify your cut criteria\n");
2012 printf("To cut on sqrt(dEta^2+dPhi^2), use: SwitchOnCutEtaPhiSum()\n");
2013 printf("To cut on dEta and dPhi separately, use: SwitchOnCutEtaPhiSeparate()\n");
2024 //------------------------------------------------------------------------------------
2025 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToEMCalSurface(AliExternalTrackParam *trkParam,
2026 const Double_t emcalR,
2027 const Double_t mass,
2028 const Double_t step,
2032 //Extrapolate track to EMCAL surface
2034 eta = -999, phi = -999;
2035 if(!trkParam) return kFALSE;
2036 if(!AliTrackerBase::PropagateTrackToBxByBz(trkParam, emcalR, mass, step, kTRUE, 0.8, -1)) return kFALSE;
2037 Double_t trkPos[3] = {0.,0.,0.};
2038 if(!trkParam->GetXYZ(trkPos)) return kFALSE;
2039 TVector3 trkPosVec(trkPos[0],trkPos[1],trkPos[2]);
2040 eta = trkPosVec.Eta();
2041 phi = trkPosVec.Phi();
2043 phi += 2*TMath::Pi();
2048 //-----------------------------------------------------------------------------------
2049 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToPosition(AliExternalTrackParam *trkParam,
2050 const Float_t *clsPos,
2057 //Return the residual by extrapolating a track param to a global position
2061 if(!trkParam) return kFALSE;
2062 Double_t trkPos[3] = {0.,0.,0.};
2063 TVector3 vec(clsPos[0],clsPos[1],clsPos[2]);
2064 Double_t alpha = ((int)(vec.Phi()*TMath::RadToDeg()/20)+0.5)*20*TMath::DegToRad();
2065 vec.RotateZ(-alpha); //Rotate the cluster to the local extrapolation coordinate system
2066 if(!AliTrackerBase::PropagateTrackToBxByBz(trkParam, vec.X(), mass, step,kTRUE, 0.8, -1)) return kFALSE;
2067 if(!trkParam->GetXYZ(trkPos)) return kFALSE; //Get the extrapolated global position
2069 TVector3 clsPosVec(clsPos[0],clsPos[1],clsPos[2]);
2070 TVector3 trkPosVec(trkPos[0],trkPos[1],trkPos[2]);
2072 // track cluster matching
2073 tmpPhi = clsPosVec.DeltaPhi(trkPosVec); // tmpPhi is between -pi and pi
2074 tmpEta = clsPosVec.Eta()-trkPosVec.Eta();
2079 //----------------------------------------------------------------------------------
2080 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToCluster(AliExternalTrackParam *trkParam,
2081 const AliVCluster *cluster,
2082 const Double_t mass,
2083 const Double_t step,
2088 //Return the residual by extrapolating a track param to a cluster
2092 if(!cluster || !trkParam) return kFALSE;
2094 Float_t clsPos[3] = {0.,0.,0.};
2095 cluster->GetPosition(clsPos);
2097 return ExtrapolateTrackToPosition(trkParam, clsPos, mass, step, tmpEta, tmpPhi);
2100 //---------------------------------------------------------------------------------
2101 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToCluster(AliExternalTrackParam *trkParam,
2102 const AliVCluster *cluster,
2107 //Return the residual by extrapolating a track param to a clusterfStepCluster
2110 return ExtrapolateTrackToCluster(trkParam, cluster, fMass, fStepCluster, tmpEta, tmpPhi);
2113 //_______________________________________________________________________
2114 void AliEMCALRecoUtils::GetMatchedResiduals(const Int_t clsIndex,
2115 Float_t &dEta, Float_t &dPhi)
2117 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2118 //Get the residuals dEta and dPhi for this cluster to the closest track
2119 //Works with ESDs and AODs
2121 if( FindMatchedPosForCluster(clsIndex) >= 999 )
2123 AliDebug(2,"No matched tracks found!\n");
2128 dEta = fResidualEta->At(FindMatchedPosForCluster(clsIndex));
2129 dPhi = fResidualPhi->At(FindMatchedPosForCluster(clsIndex));
2132 //______________________________________________________________________________________________
2133 void AliEMCALRecoUtils::GetMatchedClusterResiduals(Int_t trkIndex, Float_t &dEta, Float_t &dPhi)
2135 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2136 //Get the residuals dEta and dPhi for this track to the closest cluster
2137 //Works with ESDs and AODs
2139 if( FindMatchedPosForTrack(trkIndex) >= 999 )
2141 AliDebug(2,"No matched cluster found!\n");
2146 dEta = fResidualEta->At(FindMatchedPosForTrack(trkIndex));
2147 dPhi = fResidualPhi->At(FindMatchedPosForTrack(trkIndex));
2150 //__________________________________________________________
2151 Int_t AliEMCALRecoUtils::GetMatchedTrackIndex(Int_t clsIndex)
2153 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2154 //Get the index of matched track to this cluster
2155 //Works with ESDs and AODs
2157 if(IsClusterMatched(clsIndex))
2158 return fMatchedTrackIndex->At(FindMatchedPosForCluster(clsIndex));
2163 //__________________________________________________________
2164 Int_t AliEMCALRecoUtils::GetMatchedClusterIndex(Int_t trkIndex)
2166 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2167 //Get the index of matched cluster to this track
2168 //Works with ESDs and AODs
2170 if(IsTrackMatched(trkIndex))
2171 return fMatchedClusterIndex->At(FindMatchedPosForTrack(trkIndex));
2176 //______________________________________________________________
2177 Bool_t AliEMCALRecoUtils::IsClusterMatched(Int_t clsIndex) const
2179 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2180 //Returns if the cluster has a match
2181 if(FindMatchedPosForCluster(clsIndex) < 999)
2187 //____________________________________________________________
2188 Bool_t AliEMCALRecoUtils::IsTrackMatched(Int_t trkIndex) const
2190 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2191 //Returns if the track has a match
2192 if(FindMatchedPosForTrack(trkIndex) < 999)
2198 //______________________________________________________________________
2199 UInt_t AliEMCALRecoUtils::FindMatchedPosForCluster(Int_t clsIndex) const
2201 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2202 //Returns the position of the match in the fMatchedClusterIndex array
2203 Float_t tmpR = fCutR;
2206 for(Int_t i=0; i<fMatchedClusterIndex->GetSize(); i++)
2208 if(fMatchedClusterIndex->At(i)==clsIndex)
2210 Float_t r = TMath::Sqrt(fResidualEta->At(i)*fResidualEta->At(i) + fResidualPhi->At(i)*fResidualPhi->At(i));
2215 AliDebug(3,Form("Matched cluster index: index: %d, dEta: %2.4f, dPhi: %2.4f.\n",
2216 fMatchedClusterIndex->At(i),fResidualEta->At(i),fResidualPhi->At(i)));
2223 //____________________________________________________________________
2224 UInt_t AliEMCALRecoUtils::FindMatchedPosForTrack(Int_t trkIndex) const
2226 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2227 //Returns the position of the match in the fMatchedTrackIndex array
2228 Float_t tmpR = fCutR;
2231 for(Int_t i=0; i<fMatchedTrackIndex->GetSize(); i++)
2233 if(fMatchedTrackIndex->At(i)==trkIndex)
2235 Float_t r = TMath::Sqrt(fResidualEta->At(i)*fResidualEta->At(i) + fResidualPhi->At(i)*fResidualPhi->At(i));
2240 AliDebug(3,Form("Matched track index: index: %d, dEta: %2.4f, dPhi: %2.4f.\n",
2241 fMatchedTrackIndex->At(i),fResidualEta->At(i),fResidualPhi->At(i)));
2248 //__________________________________________________________________________
2249 Bool_t AliEMCALRecoUtils::IsGoodCluster(AliVCluster *cluster,
2250 const AliEMCALGeometry *geom,
2251 AliVCaloCells* cells,const Int_t bc)
2253 // check if the cluster survives some quality cut
2256 Bool_t isGood=kTRUE;
2258 if(!cluster || !cluster->IsEMCAL()) return kFALSE;
2260 if(ClusterContainsBadChannel(geom,cluster->GetCellsAbsId(),cluster->GetNCells())) return kFALSE;
2262 if(!CheckCellFiducialRegion(geom,cluster,cells)) return kFALSE;
2264 if(IsExoticCluster(cluster, cells,bc)) return kFALSE;
2269 //__________________________________________________________
2270 Bool_t AliEMCALRecoUtils::IsAccepted(AliESDtrack *esdTrack)
2272 // Given a esd track, return whether the track survive all the cuts
2274 // The different quality parameter are first
2275 // retrieved from the track. then it is found out what cuts the
2276 // track did not survive and finally the cuts are imposed.
2278 UInt_t status = esdTrack->GetStatus();
2280 Int_t nClustersITS = esdTrack->GetITSclusters(0);
2281 Int_t nClustersTPC = esdTrack->GetTPCclusters(0);
2283 Float_t chi2PerClusterITS = -1;
2284 Float_t chi2PerClusterTPC = -1;
2285 if (nClustersITS!=0)
2286 chi2PerClusterITS = esdTrack->GetITSchi2()/Float_t(nClustersITS);
2287 if (nClustersTPC!=0)
2288 chi2PerClusterTPC = esdTrack->GetTPCchi2()/Float_t(nClustersTPC);
2292 if(fTrackCutsType==kGlobalCut)
2294 Float_t maxDCAToVertexXYPtDep = 0.0182 + 0.0350/TMath::Power(esdTrack->Pt(),1.01); //This expression comes from AliESDtrackCuts::GetStandardITSTPCTrackCuts2010()
2295 //AliDebug(3,Form("Track pT = %f, DCAtoVertexXY = %f",esdTrack->Pt(),MaxDCAToVertexXYPtDep));
2296 SetMaxDCAToVertexXY(maxDCAToVertexXYPtDep); //Set pT dependent DCA cut to vertex in x-y plane
2302 esdTrack->GetImpactParameters(b,bCov);
2303 if (bCov[0]<=0 || bCov[2]<=0)
2305 AliDebug(1, "Estimated b resolution lower or equal zero!");
2306 bCov[0]=0; bCov[2]=0;
2309 Float_t dcaToVertexXY = b[0];
2310 Float_t dcaToVertexZ = b[1];
2311 Float_t dcaToVertex = -1;
2313 if (fCutDCAToVertex2D)
2314 dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY/fCutMaxDCAToVertexXY/fCutMaxDCAToVertexXY + dcaToVertexZ*dcaToVertexZ/fCutMaxDCAToVertexZ/fCutMaxDCAToVertexZ);
2316 dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY + dcaToVertexZ*dcaToVertexZ);
2320 Bool_t cuts[kNCuts];
2321 for (Int_t i=0; i<kNCuts; i++) cuts[i]=kFALSE;
2323 // track quality cuts
2324 if (fCutRequireTPCRefit && (status&AliESDtrack::kTPCrefit)==0)
2326 if (fCutRequireITSRefit && (status&AliESDtrack::kITSrefit)==0)
2328 if (nClustersTPC<fCutMinNClusterTPC)
2330 if (nClustersITS<fCutMinNClusterITS)
2332 if (chi2PerClusterTPC>fCutMaxChi2PerClusterTPC)
2334 if (chi2PerClusterITS>fCutMaxChi2PerClusterITS)
2336 if (!fCutAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0)
2338 if (fCutDCAToVertex2D && dcaToVertex > 1)
2340 if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexXY) > fCutMaxDCAToVertexXY)
2342 if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexZ) > fCutMaxDCAToVertexZ)
2345 if(fTrackCutsType==kGlobalCut)
2347 //Require at least one SPD point + anything else in ITS
2348 if( (esdTrack->HasPointOnITSLayer(0) || esdTrack->HasPointOnITSLayer(1)) == kFALSE)
2353 if(fCutRequireITSStandAlone || fCutRequireITSpureSA){
2354 if ((status & AliESDtrack::kITSin) == 0 || (status & AliESDtrack::kTPCin)){
2358 // ITS standalone tracks
2359 if(fCutRequireITSStandAlone && !fCutRequireITSpureSA){
2360 if(status & AliESDtrack::kITSpureSA) cuts[11] = kTRUE;
2361 }else if(fCutRequireITSpureSA){
2362 if(!(status & AliESDtrack::kITSpureSA)) cuts[11] = kTRUE;
2368 for (Int_t i=0; i<kNCuts; i++)
2369 if (cuts[i]) { cut = kTRUE ; }
2378 //_____________________________________
2379 void AliEMCALRecoUtils::InitTrackCuts()
2381 //Intilize the track cut criteria
2382 //By default these cuts are set according to AliESDtrackCuts::GetStandardTPCOnlyTrackCuts()
2383 //Also you can customize the cuts using the setters
2385 switch (fTrackCutsType)
2389 AliInfo(Form("Track cuts for matching: GetStandardTPCOnlyTrackCuts()"));
2391 SetMinNClustersTPC(70);
2392 SetMaxChi2PerClusterTPC(4);
2393 SetAcceptKinkDaughters(kFALSE);
2394 SetRequireTPCRefit(kFALSE);
2397 SetRequireITSRefit(kFALSE);
2398 SetMaxDCAToVertexZ(3.2);
2399 SetMaxDCAToVertexXY(2.4);
2400 SetDCAToVertex2D(kTRUE);
2407 AliInfo(Form("Track cuts for matching: GetStandardITSTPCTrackCuts2010(kTURE)"));
2409 SetMinNClustersTPC(70);
2410 SetMaxChi2PerClusterTPC(4);
2411 SetAcceptKinkDaughters(kFALSE);
2412 SetRequireTPCRefit(kTRUE);
2415 SetRequireITSRefit(kTRUE);
2416 SetMaxDCAToVertexZ(2);
2417 SetMaxDCAToVertexXY();
2418 SetDCAToVertex2D(kFALSE);
2425 AliInfo(Form("Track cuts for matching: Loose cut w/o DCA cut"));
2426 SetMinNClustersTPC(50);
2427 SetAcceptKinkDaughters(kTRUE);
2432 case kITSStandAlone:
2434 AliInfo(Form("Track cuts for matching: ITS Stand Alone tracks cut w/o DCA cut"));
2435 SetRequireITSRefit(kTRUE);
2436 SetRequireITSStandAlone(kTRUE);
2437 SetITSTrackSA(kTRUE);
2445 //________________________________________________________________________
2446 void AliEMCALRecoUtils::SetClusterMatchedToTrack(const AliVEvent *event)
2448 // Checks if tracks are matched to EMC clusters and set the matched EMCAL cluster index to ESD track.
2450 Int_t nTracks = event->GetNumberOfTracks();
2451 for (Int_t iTrack = 0; iTrack < nTracks; ++iTrack)
2453 AliVTrack* track = dynamic_cast<AliVTrack*>(event->GetTrack(iTrack));
2456 AliWarning(Form("Could not receive track %d", iTrack));
2460 Int_t matchClusIndex = GetMatchedClusterIndex(iTrack);
2461 track->SetEMCALcluster(matchClusIndex); //sets -1 if track not matched within residual
2462 /*the following can be done better if AliVTrack::SetStatus will be there. Patch pending with Andreas/Peter*/
2463 AliESDtrack* esdtrack = dynamic_cast<AliESDtrack*>(track);
2465 if(matchClusIndex != -1)
2466 esdtrack->SetStatus(AliESDtrack::kEMCALmatch);
2468 esdtrack->ResetStatus(AliESDtrack::kEMCALmatch);
2470 AliAODTrack* aodtrack = dynamic_cast<AliAODTrack*>(track);
2471 if(matchClusIndex != -1)
2472 aodtrack->SetStatus(AliESDtrack::kEMCALmatch);
2474 aodtrack->ResetStatus(AliESDtrack::kEMCALmatch);
2478 AliDebug(2,"Track matched to closest cluster");
2481 //_________________________________________________________________________
2482 void AliEMCALRecoUtils::SetTracksMatchedToCluster(const AliVEvent *event)
2484 // Checks if EMC clusters are matched to ESD track.
2485 // Adds track indexes of all the tracks matched to a cluster withing residuals in ESDCalocluster.
2487 for (Int_t iClus=0; iClus < event->GetNumberOfCaloClusters(); ++iClus)
2489 AliVCluster *cluster = event->GetCaloCluster(iClus);
2490 if (!cluster->IsEMCAL())
2493 Int_t nTracks = event->GetNumberOfTracks();
2494 TArrayI arrayTrackMatched(nTracks);
2496 // Get the closest track matched to the cluster
2498 Int_t matchTrackIndex = GetMatchedTrackIndex(iClus);
2499 if (matchTrackIndex != -1)
2501 arrayTrackMatched[nMatched] = matchTrackIndex;
2505 // Get all other tracks matched to the cluster
2506 for(Int_t iTrk=0; iTrk<nTracks; ++iTrk)
2508 AliVTrack* track = dynamic_cast<AliVTrack*>(event->GetTrack(iTrk));
2509 if(iTrk == matchTrackIndex) continue;
2510 if(track->GetEMCALcluster() == iClus)
2512 arrayTrackMatched[nMatched] = iTrk;
2517 //printf("Tender::SetTracksMatchedToCluster - cluster E %f, N matches %d, first match %d\n",cluster->E(),nMatched,arrayTrackMatched[0]);
2519 arrayTrackMatched.Set(nMatched);
2520 AliESDCaloCluster *esdcluster = dynamic_cast<AliESDCaloCluster*>(cluster);
2522 esdcluster->AddTracksMatched(arrayTrackMatched);
2523 else if (nMatched>0) {
2524 AliAODCaloCluster *aodcluster = dynamic_cast<AliAODCaloCluster*>(cluster);
2526 aodcluster->AddTrackMatched(event->GetTrack(arrayTrackMatched.At(0)));
2529 Float_t eta= -999, phi = -999;
2530 if (matchTrackIndex != -1)
2531 GetMatchedResiduals(iClus, eta, phi);
2532 cluster->SetTrackDistance(phi, eta);
2535 AliDebug(2,"Cluster matched to tracks");
2538 //___________________________________________________
2539 void AliEMCALRecoUtils::Print(const Option_t *) const
2543 printf("AliEMCALRecoUtils Settings: \n");
2544 printf("Misalignment shifts\n");
2545 for(Int_t i=0; i<5; i++) printf("\t sector %d, traslation (x,y,z)=(%f,%f,%f), rotation (x,y,z)=(%f,%f,%f)\n",i,
2546 fMisalTransShift[i*3],fMisalTransShift[i*3+1],fMisalTransShift[i*3+2],
2547 fMisalRotShift[i*3], fMisalRotShift[i*3+1], fMisalRotShift[i*3+2] );
2548 printf("Non linearity function %d, parameters:\n", fNonLinearityFunction);
2549 for(Int_t i=0; i<6; i++) printf("param[%d]=%f\n",i, fNonLinearityParams[i]);
2551 printf("Position Recalculation option %d, Particle Type %d, fW0 %2.2f, Recalibrate Data %d \n",fPosAlgo,fParticleType,fW0, fRecalibration);
2553 printf("Matching criteria: ");
2556 printf("sqrt(dEta^2+dPhi^2)<%4.3f\n",fCutR);
2558 else if(fCutEtaPhiSeparate)
2560 printf("dEta<%4.3f, dPhi<%4.3f\n",fCutEta,fCutPhi);
2565 printf("please specify your cut criteria\n");
2566 printf("To cut on sqrt(dEta^2+dPhi^2), use: SwitchOnCutEtaPhiSum()\n");
2567 printf("To cut on dEta and dPhi separately, use: SwitchOnCutEtaPhiSeparate()\n");
2570 printf("Mass hypothesis = %2.3f [GeV/c^2], extrapolation step to surface = %2.2f[cm], step to cluster = %2.2f[cm]\n",fMass,fStepSurface, fStepCluster);
2571 printf("Cluster selection window: dR < %2.0f\n",fClusterWindow);
2573 printf("Track cuts: \n");
2574 printf("Minimum track pT: %1.2f\n",fCutMinTrackPt);
2575 printf("AOD track selection mask: %d\n",fAODFilterMask);
2576 printf("TPCRefit = %d, ITSRefit = %d\n",fCutRequireTPCRefit,fCutRequireITSRefit);
2577 printf("AcceptKinks = %d\n",fCutAcceptKinkDaughters);
2578 printf("MinNCulsterTPC = %d, MinNClusterITS = %d\n",fCutMinNClusterTPC,fCutMinNClusterITS);
2579 printf("MaxChi2TPC = %2.2f, MaxChi2ITS = %2.2f\n",fCutMaxChi2PerClusterTPC,fCutMaxChi2PerClusterITS);
2580 printf("DCSToVertex2D = %d, MaxDCAToVertexXY = %2.2f, MaxDCAToVertexZ = %2.2f\n",fCutDCAToVertex2D,fCutMaxDCAToVertexXY,fCutMaxDCAToVertexZ);