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 | Sun Dec 8 06:56:48 2013 +0100 | Constantin Loizides $ */
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 ClassImp(AliEMCALRecoUtils)
63 //_____________________________________
64 AliEMCALRecoUtils::AliEMCALRecoUtils():
65 fParticleType(0), fPosAlgo(0), fW0(0),
66 fNonLinearityFunction(0), fNonLinearThreshold(0),
67 fSmearClusterEnergy(kFALSE), fRandom(),
68 fCellsRecalibrated(kFALSE), fRecalibration(kFALSE), fEMCALRecalibrationFactors(),
69 fTimeRecalibration(kFALSE), fEMCALTimeRecalibrationFactors(), fUseRunCorrectionFactors(kFALSE),
70 fRemoveBadChannels(kFALSE), fRecalDistToBadChannels(kFALSE), fEMCALBadChannelMap(),
71 fNCellsFromEMCALBorder(0), fNoEMCALBorderAtEta0(kTRUE),
72 fRejectExoticCluster(kFALSE), fRejectExoticCells(kFALSE),
73 fExoticCellFraction(0), fExoticCellDiffTime(0), fExoticCellMinAmplitude(0),
74 fPIDUtils(), fAODFilterMask(0),
75 fAODHybridTracks(0), fAODTPCOnlyTracks(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),
81 fITSTrackSA(kFALSE), fEMCalSurfaceDistance(440.),
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 fAODHybridTracks(reco.fAODHybridTracks), fAODTPCOnlyTracks(reco.fAODTPCOnlyTracks),
124 fMatchedTrackIndex( reco.fMatchedTrackIndex? new TArrayI(*reco.fMatchedTrackIndex):0x0),
125 fMatchedClusterIndex(reco.fMatchedClusterIndex?new TArrayI(*reco.fMatchedClusterIndex):0x0),
126 fResidualEta( reco.fResidualEta? new TArrayF(*reco.fResidualEta):0x0),
127 fResidualPhi( reco.fResidualPhi? new TArrayF(*reco.fResidualPhi):0x0),
128 fCutEtaPhiSum(reco.fCutEtaPhiSum), fCutEtaPhiSeparate(reco.fCutEtaPhiSeparate),
129 fCutR(reco.fCutR), fCutEta(reco.fCutEta), fCutPhi(reco.fCutPhi),
130 fClusterWindow(reco.fClusterWindow),
131 fMass(reco.fMass), fStepSurface(reco.fStepSurface), fStepCluster(reco.fStepCluster),
132 fITSTrackSA(reco.fITSTrackSA), fEMCalSurfaceDistance(440.),
133 fTrackCutsType(reco.fTrackCutsType), fCutMinTrackPt(reco.fCutMinTrackPt),
134 fCutMinNClusterTPC(reco.fCutMinNClusterTPC), fCutMinNClusterITS(reco.fCutMinNClusterITS),
135 fCutMaxChi2PerClusterTPC(reco.fCutMaxChi2PerClusterTPC), fCutMaxChi2PerClusterITS(reco.fCutMaxChi2PerClusterITS),
136 fCutRequireTPCRefit(reco.fCutRequireTPCRefit), fCutRequireITSRefit(reco.fCutRequireITSRefit),
137 fCutAcceptKinkDaughters(reco.fCutAcceptKinkDaughters), fCutMaxDCAToVertexXY(reco.fCutMaxDCAToVertexXY),
138 fCutMaxDCAToVertexZ(reco.fCutMaxDCAToVertexZ), fCutDCAToVertex2D(reco.fCutDCAToVertex2D),
139 fCutRequireITSStandAlone(reco.fCutRequireITSStandAlone), fCutRequireITSpureSA(reco.fCutRequireITSpureSA)
143 for (Int_t i = 0; i < 15 ; i++) { fMisalRotShift[i] = reco.fMisalRotShift[i] ;
144 fMisalTransShift[i] = reco.fMisalTransShift[i] ; }
145 for (Int_t i = 0; i < 7 ; i++) { fNonLinearityParams[i] = reco.fNonLinearityParams[i] ; }
146 for (Int_t i = 0; i < 3 ; i++) { fSmearClusterParam[i] = reco.fSmearClusterParam[i] ; }
148 fUseMassForTracking = reco.fUseMassForTracking;
153 //______________________________________________________________________
154 AliEMCALRecoUtils & AliEMCALRecoUtils::operator = (const AliEMCALRecoUtils & reco)
156 //Assignment operator
158 if (this == &reco)return *this;
159 ((TNamed *)this)->operator=(reco);
161 for (Int_t i = 0; i < 15 ; i++) { fMisalTransShift[i] = reco.fMisalTransShift[i] ;
162 fMisalRotShift[i] = reco.fMisalRotShift[i] ; }
163 for (Int_t i = 0; i < 7 ; i++) { fNonLinearityParams[i] = reco.fNonLinearityParams[i] ; }
164 for (Int_t i = 0; i < 3 ; i++) { fSmearClusterParam[i] = reco.fSmearClusterParam[i] ; }
166 fParticleType = reco.fParticleType;
167 fPosAlgo = reco.fPosAlgo;
170 fNonLinearityFunction = reco.fNonLinearityFunction;
171 fNonLinearThreshold = reco.fNonLinearThreshold;
172 fSmearClusterEnergy = reco.fSmearClusterEnergy;
174 fCellsRecalibrated = reco.fCellsRecalibrated;
175 fRecalibration = reco.fRecalibration;
176 fEMCALRecalibrationFactors = reco.fEMCALRecalibrationFactors;
178 fTimeRecalibration = reco.fTimeRecalibration;
179 fEMCALTimeRecalibrationFactors = reco.fEMCALTimeRecalibrationFactors;
181 fUseRunCorrectionFactors = reco.fUseRunCorrectionFactors;
183 fRemoveBadChannels = reco.fRemoveBadChannels;
184 fRecalDistToBadChannels = reco.fRecalDistToBadChannels;
185 fEMCALBadChannelMap = reco.fEMCALBadChannelMap;
187 fNCellsFromEMCALBorder = reco.fNCellsFromEMCALBorder;
188 fNoEMCALBorderAtEta0 = reco.fNoEMCALBorderAtEta0;
190 fRejectExoticCluster = reco.fRejectExoticCluster;
191 fRejectExoticCells = reco.fRejectExoticCells;
192 fExoticCellFraction = reco.fExoticCellFraction;
193 fExoticCellDiffTime = reco.fExoticCellDiffTime;
194 fExoticCellMinAmplitude = reco.fExoticCellMinAmplitude;
196 fPIDUtils = reco.fPIDUtils;
198 fAODFilterMask = reco.fAODFilterMask;
199 fAODHybridTracks = reco.fAODHybridTracks;
200 fAODTPCOnlyTracks = reco.fAODTPCOnlyTracks;
202 fCutEtaPhiSum = reco.fCutEtaPhiSum;
203 fCutEtaPhiSeparate = reco.fCutEtaPhiSeparate;
205 fCutEta = reco.fCutEta;
206 fCutPhi = reco.fCutPhi;
207 fClusterWindow = reco.fClusterWindow;
209 fStepSurface = reco.fStepSurface;
210 fStepCluster = reco.fStepCluster;
211 fITSTrackSA = reco.fITSTrackSA;
212 fEMCalSurfaceDistance = reco.fEMCalSurfaceDistance;
214 fTrackCutsType = reco.fTrackCutsType;
215 fCutMinTrackPt = reco.fCutMinTrackPt;
216 fCutMinNClusterTPC = reco.fCutMinNClusterTPC;
217 fCutMinNClusterITS = reco.fCutMinNClusterITS;
218 fCutMaxChi2PerClusterTPC = reco.fCutMaxChi2PerClusterTPC;
219 fCutMaxChi2PerClusterITS = reco.fCutMaxChi2PerClusterITS;
220 fCutRequireTPCRefit = reco.fCutRequireTPCRefit;
221 fCutRequireITSRefit = reco.fCutRequireITSRefit;
222 fCutAcceptKinkDaughters = reco.fCutAcceptKinkDaughters;
223 fCutMaxDCAToVertexXY = reco.fCutMaxDCAToVertexXY;
224 fCutMaxDCAToVertexZ = reco.fCutMaxDCAToVertexZ;
225 fCutDCAToVertex2D = reco.fCutDCAToVertex2D;
226 fCutRequireITSStandAlone = reco.fCutRequireITSStandAlone;
227 fCutRequireITSpureSA = reco.fCutRequireITSpureSA;
228 fUseMassForTracking = reco.fUseMassForTracking;
230 if (reco.fResidualEta) {
231 // assign or copy construct
233 *fResidualEta = *reco.fResidualEta;
235 fResidualEta = new TArrayF(*reco.fResidualEta);
238 if (fResidualEta) delete fResidualEta;
242 if (reco.fResidualPhi) {
243 // assign or copy construct
245 *fResidualPhi = *reco.fResidualPhi;
247 fResidualPhi = new TArrayF(*reco.fResidualPhi);
250 if (fResidualPhi) delete fResidualPhi;
254 if (reco.fMatchedTrackIndex) {
255 // assign or copy construct
256 if (fMatchedTrackIndex) {
257 *fMatchedTrackIndex = *reco.fMatchedTrackIndex;
259 fMatchedTrackIndex = new TArrayI(*reco.fMatchedTrackIndex);
262 if (fMatchedTrackIndex) delete fMatchedTrackIndex;
263 fMatchedTrackIndex = 0;
266 if (reco.fMatchedClusterIndex) {
267 // assign or copy construct
268 if (fMatchedClusterIndex) {
269 *fMatchedClusterIndex = *reco.fMatchedClusterIndex;
271 fMatchedClusterIndex = new TArrayI(*reco.fMatchedClusterIndex);
274 if (fMatchedClusterIndex) delete fMatchedClusterIndex;
275 fMatchedClusterIndex = 0;
281 //_____________________________________
282 AliEMCALRecoUtils::~AliEMCALRecoUtils()
286 if (fEMCALRecalibrationFactors) {
287 fEMCALRecalibrationFactors->Clear();
288 delete fEMCALRecalibrationFactors;
291 if (fEMCALTimeRecalibrationFactors) {
292 fEMCALTimeRecalibrationFactors->Clear();
293 delete fEMCALTimeRecalibrationFactors;
296 if (fEMCALBadChannelMap) {
297 fEMCALBadChannelMap->Clear();
298 delete fEMCALBadChannelMap;
301 delete fMatchedTrackIndex ;
302 delete fMatchedClusterIndex ;
303 delete fResidualEta ;
304 delete fResidualPhi ;
310 //_______________________________________________________________________________
311 Bool_t AliEMCALRecoUtils::AcceptCalibrateCell(Int_t absID, Int_t bc,
312 Float_t & amp, Double_t & time,
313 AliVCaloCells* cells)
315 // Reject cell if criteria not passed and calibrate it
317 AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
319 if (absID < 0 || absID >= 24*48*geom->GetNumberOfSuperModules())
322 Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1;
324 if (!geom->GetCellIndex(absID,imod,iTower,iIphi,iIeta)) {
325 // cell absID does not exist
330 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);
332 // Do not include bad channels found in analysis,
333 if (IsBadChannelsRemovalSwitchedOn() && GetEMCALChannelStatus(imod, ieta, iphi)) {
338 amp = cells->GetCellAmplitude(absID);
339 if (!fCellsRecalibrated && IsRecalibrationOn())
340 amp *= GetEMCALChannelRecalibrationFactor(imod,ieta,iphi);
343 time = cells->GetCellTime(absID);
345 RecalibrateCellTime(absID,bc,time);
350 //_____________________________________________________________________________
351 Bool_t AliEMCALRecoUtils::CheckCellFiducialRegion(const AliEMCALGeometry* geom,
352 const AliVCluster* cluster,
353 AliVCaloCells* cells)
355 // Given the list of AbsId of the cluster, get the maximum cell and
356 // check if there are fNCellsFromBorder from the calorimeter border
360 AliInfo("Cluster pointer null!");
364 //If the distance to the border is 0 or negative just exit accept all clusters
365 if (cells->GetType()==AliVCaloCells::kEMCALCell && fNCellsFromEMCALBorder <= 0 )
368 Int_t absIdMax = -1, iSM =-1, ieta = -1, iphi = -1;
369 Bool_t shared = kFALSE;
370 GetMaxEnergyCell(geom, cells, cluster, absIdMax, iSM, ieta, iphi, shared);
372 AliDebug(2,Form("Cluster Max AbsId %d, Cell Energy %2.2f, Cluster Energy %2.2f, Ncells from border %d, EMCAL eta=0 %d\n",
373 absIdMax, cells->GetCellAmplitude(absIdMax), cluster->E(), fNCellsFromEMCALBorder, fNoEMCALBorderAtEta0));
375 if (absIdMax==-1) return kFALSE;
377 //Check if the cell is close to the borders:
378 Bool_t okrow = kFALSE;
379 Bool_t okcol = kFALSE;
381 if (iSM < 0 || iphi < 0 || ieta < 0 ) {
382 AliFatal(Form("Negative value for super module: %d, or cell ieta: %d, or cell iphi: %d, check EMCAL geometry name\n",
384 return kFALSE; // trick coverity
389 if( geom->GetSMType(iSM) == AliEMCALGeometry::kEMCAL_Half ) iPhiLast /= 2;
390 else if ( geom->GetSMType(iSM) == AliEMCALGeometry::kEMCAL_3rd ) iPhiLast /= 3;// 1/3 sm case
392 if(iphi >= fNCellsFromEMCALBorder && iphi < iPhiLast - fNCellsFromEMCALBorder) okrow = kTRUE;
396 if(!fNoEMCALBorderAtEta0 || geom->IsDCALSM(iSM)) {// conside inner border
397 if( geom->GetSMType(iSM) == AliEMCALGeometry::kDCAL_Standard ) iEtaLast = iEtaLast*2/3;
398 if(ieta > fNCellsFromEMCALBorder && ieta < iEtaLast-fNCellsFromEMCALBorder) okcol = kTRUE;
401 if (ieta >= fNCellsFromEMCALBorder) okcol = kTRUE;
403 if(ieta < iEtaLast-fNCellsFromEMCALBorder) okcol = kTRUE;
407 AliDebug(2,Form("EMCAL Cluster in %d cells fiducial volume: ieta %d, iphi %d, SM %d: column? %d, row? %d\nq",
408 fNCellsFromEMCALBorder, ieta, iphi, iSM, okcol, okrow));
410 if (okcol && okrow) {
411 //printf("Accept\n");
414 //printf("Reject\n");
415 AliDebug(2,Form("Reject cluster in border, max cell : ieta %d, iphi %d, SM %d\n",ieta, iphi, iSM));
420 //_______________________________________________________________________________
421 Bool_t AliEMCALRecoUtils::ClusterContainsBadChannel(const AliEMCALGeometry* geom,
422 const UShort_t* cellList,
425 // Check that in the cluster cells, there is no bad channel of those stored
426 // in fEMCALBadChannelMap or fPHOSBadChannelMap
428 if (!fRemoveBadChannels) return kFALSE;
429 if (!fEMCALBadChannelMap) return kFALSE;
434 for (Int_t iCell = 0; iCell<nCells; iCell++) {
435 //Get the column and row
436 Int_t iTower = -1, iIphi = -1, iIeta = -1;
437 geom->GetCellIndex(cellList[iCell],imod,iTower,iIphi,iIeta);
438 if (fEMCALBadChannelMap->GetEntries() <= imod) continue;
439 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,irow,icol);
440 if (GetEMCALChannelStatus(imod, icol, irow)) {
441 AliDebug(2,Form("Cluster with bad channel: SM %d, col %d, row %d\n",imod, icol, irow));
444 }// cell cluster loop
450 //___________________________________________________________________________
451 Float_t AliEMCALRecoUtils::GetECross(Int_t absID, Double_t tcell,
452 AliVCaloCells* cells, Int_t bc)
454 //Calculate the energy in the cross around the energy given cell
456 AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
458 Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1;
459 geom->GetCellIndex(absID,imod,iTower,iIphi,iIeta);
460 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);
462 //Get close cells index, energy and time, not in corners
467 if ( iphi < AliEMCALGeoParams::fgkEMCALRows-1) absID1 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi+1, ieta);
468 if ( iphi > 0 ) absID2 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi-1, ieta);
470 // In case of cell in eta = 0 border, depending on SM shift the cross cell index
475 if ( ieta == AliEMCALGeoParams::fgkEMCALCols-1 && !(imod%2) ) {
476 absID3 = geom-> GetAbsCellIdFromCellIndexes(imod+1, iphi, 0);
477 absID4 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta-1);
478 } else if ( ieta == 0 && imod%2 ) {
479 absID3 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta+1);
480 absID4 = geom-> GetAbsCellIdFromCellIndexes(imod-1, iphi, AliEMCALGeoParams::fgkEMCALCols-1);
482 if ( ieta < AliEMCALGeoParams::fgkEMCALCols-1 )
483 absID3 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta+1);
485 absID4 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta-1);
488 //printf("IMOD %d, AbsId %d, a %d, b %d, c %d e %d \n",imod,absID,absID1,absID2,absID3,absID4);
490 Float_t ecell1 = 0, ecell2 = 0, ecell3 = 0, ecell4 = 0;
491 Double_t tcell1 = 0, tcell2 = 0, tcell3 = 0, tcell4 = 0;
493 AcceptCalibrateCell(absID1,bc, ecell1,tcell1,cells);
494 AcceptCalibrateCell(absID2,bc, ecell2,tcell2,cells);
495 AcceptCalibrateCell(absID3,bc, ecell3,tcell3,cells);
496 AcceptCalibrateCell(absID4,bc, ecell4,tcell4,cells);
498 if (TMath::Abs(tcell-tcell1)*1.e9 > fExoticCellDiffTime) ecell1 = 0 ;
499 if (TMath::Abs(tcell-tcell2)*1.e9 > fExoticCellDiffTime) ecell2 = 0 ;
500 if (TMath::Abs(tcell-tcell3)*1.e9 > fExoticCellDiffTime) ecell3 = 0 ;
501 if (TMath::Abs(tcell-tcell4)*1.e9 > fExoticCellDiffTime) ecell4 = 0 ;
503 return ecell1+ecell2+ecell3+ecell4;
506 //_____________________________________________________________________________________________
507 Bool_t AliEMCALRecoUtils::IsExoticCell(Int_t absID, AliVCaloCells* cells, Int_t bc)
509 // Look to cell neighbourhood and reject if it seems exotic
510 // Do before recalibrating the cells
512 if (!fRejectExoticCells) return kFALSE;
516 Bool_t accept = AcceptCalibrateCell(absID, bc, ecell ,tcell ,cells);
518 if (!accept) return kTRUE; // reject this cell
520 if (ecell < fExoticCellMinAmplitude) return kFALSE; // do not reject low energy cells
522 Float_t eCross = GetECross(absID,tcell,cells,bc);
524 if (1-eCross/ecell > fExoticCellFraction) {
525 AliDebug(2,Form("AliEMCALRecoUtils::IsExoticCell() - EXOTIC CELL id %d, eCell %f, eCross %f, 1-eCross/eCell %f\n",
526 absID,ecell,eCross,1-eCross/ecell));
533 //___________________________________________________________________
534 Bool_t AliEMCALRecoUtils::IsExoticCluster(const AliVCluster *cluster,
535 AliVCaloCells *cells,
538 // Check if the cluster highest energy tower is exotic
541 AliInfo("Cluster pointer null!");
545 if (!fRejectExoticCluster) return kFALSE;
547 // Get highest energy tower
548 AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
549 Int_t iSupMod = -1, absId = -1, ieta = -1, iphi = -1;
550 Bool_t shared = kFALSE;
551 GetMaxEnergyCell(geom, cells, cluster, absId, iSupMod, ieta, iphi, shared);
553 return IsExoticCell(absId,cells,bc);
556 //_______________________________________________________________________
557 Float_t AliEMCALRecoUtils::SmearClusterEnergy(const AliVCluster* cluster)
559 //In case of MC analysis, smear energy to match resolution/calibration in real data
562 AliInfo("Cluster pointer null!");
566 Float_t energy = cluster->E() ;
567 Float_t rdmEnergy = energy ;
568 if (fSmearClusterEnergy) {
569 rdmEnergy = fRandom.Gaus(energy,fSmearClusterParam[0] * TMath::Sqrt(energy) +
570 fSmearClusterParam[1] * energy +
571 fSmearClusterParam[2] );
572 AliDebug(2, Form("Energy: original %f, smeared %f\n", energy, rdmEnergy));
578 //____________________________________________________________________________
579 Float_t AliEMCALRecoUtils::CorrectClusterEnergyLinearity(AliVCluster* cluster)
581 // Correct cluster energy from non linearity functions
584 AliInfo("Cluster pointer null!");
588 Float_t energy = cluster->E();
591 // Clusters with less than 50 MeV or negative are not possible
592 AliInfo(Form("Too Low Cluster energy!, E = %f < 0.05 GeV",energy));
596 switch (fNonLinearityFunction)
600 //Non-Linearity correction (from MC with function ([0]*exp(-[1]/E))+(([2]/([3]*2.*TMath::Pi())*exp(-(E-[4])^2/(2.*[3]^2)))))
601 //fNonLinearityParams[0] = 1.014;
602 //fNonLinearityParams[1] =-0.03329;
603 //fNonLinearityParams[2] =-0.3853;
604 //fNonLinearityParams[3] = 0.5423;
605 //fNonLinearityParams[4] =-0.4335;
606 energy *= (fNonLinearityParams[0]*exp(-fNonLinearityParams[1]/energy))+
607 ((fNonLinearityParams[2]/(fNonLinearityParams[3]*2.*TMath::Pi())*
608 exp(-(energy-fNonLinearityParams[4])*(energy-fNonLinearityParams[4])/(2.*fNonLinearityParams[3]*fNonLinearityParams[3]))));
614 //Non-Linearity correction (from MC with function [0]/((x+[1])^[2]))+1;
615 //fNonLinearityParams[0] = 3.11111e-02;
616 //fNonLinearityParams[1] =-5.71666e-02;
617 //fNonLinearityParams[2] = 5.67995e-01;
619 energy *= fNonLinearityParams[0]/TMath::Power(energy+fNonLinearityParams[1],fNonLinearityParams[2])+1;
625 //Same as beam test corrected, change parameters
626 //fNonLinearityParams[0] = 9.81039e-01
627 //fNonLinearityParams[1] = 1.13508e-01;
628 //fNonLinearityParams[2] = 1.00173e+00;
629 //fNonLinearityParams[3] = 9.67998e-02;
630 //fNonLinearityParams[4] = 2.19381e+02;
631 //fNonLinearityParams[5] = 6.31604e+01;
632 //fNonLinearityParams[6] = 1;
633 energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
641 //Non-Linearity correction (from Olga Data with function p0+p1*exp(-p2*E))
642 //fNonLinearityParams[0] = 1.04;
643 //fNonLinearityParams[1] = -0.1445;
644 //fNonLinearityParams[2] = 1.046;
645 energy /= (fNonLinearityParams[0]+fNonLinearityParams[1]*exp(-fNonLinearityParams[2]*energy)); //Olga function
649 case kPi0GammaConversion:
651 //Non-Linearity correction (Nicolas from Dimitri Data with function C*[1-a*exp(-b*E)])
652 //fNonLinearityParams[0] = 0.139393/0.1349766;
653 //fNonLinearityParams[1] = 0.0566186;
654 //fNonLinearityParams[2] = 0.982133;
655 energy /= fNonLinearityParams[0]*(1-fNonLinearityParams[1]*exp(-fNonLinearityParams[2]*energy));
662 //From beam test, Alexei's results, for different ZS thresholds
663 // th=30 MeV; th = 45 MeV; th = 75 MeV
664 //fNonLinearityParams[0] = 1.007; 1.003; 1.002
665 //fNonLinearityParams[1] = 0.894; 0.719; 0.797
666 //fNonLinearityParams[2] = 0.246; 0.334; 0.358
667 //Rescale the param[0] with 1.03
668 energy /= fNonLinearityParams[0]/(1+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]));
673 case kBeamTestCorrected:
675 //From beam test, corrected for material between beam and EMCAL
676 //fNonLinearityParams[0] = 0.99078
677 //fNonLinearityParams[1] = 0.161499;
678 //fNonLinearityParams[2] = 0.655166;
679 //fNonLinearityParams[3] = 0.134101;
680 //fNonLinearityParams[4] = 163.282;
681 //fNonLinearityParams[5] = 23.6904;
682 //fNonLinearityParams[6] = 0.978;
683 energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
688 case kBeamTestCorrectedv2:
690 //From beam test, corrected for material between beam and EMCAL
691 //fNonLinearityParams[0] = 0.983504;
692 //fNonLinearityParams[1] = 0.210106;
693 //fNonLinearityParams[2] = 0.897274;
694 //fNonLinearityParams[3] = 0.0829064;
695 //fNonLinearityParams[4] = 152.299;
696 //fNonLinearityParams[5] = 31.5028;
697 //fNonLinearityParams[6] = 0.968;
698 energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
705 //Based on fit to the MC/data using kNoCorrection on the data - utilizes symmetric decay method and kPi0MCv5(MC) - 28 Oct 2013
706 //fNonLinearityParams[0] = 1.0;
707 //fNonLinearityParams[1] = 6.64778e-02;
708 //fNonLinearityParams[2] = 1.570;
709 //fNonLinearityParams[3] = 9.67998e-02;
710 //fNonLinearityParams[4] = 2.19381e+02;
711 //fNonLinearityParams[5] = 6.31604e+01;
712 //fNonLinearityParams[6] = 1.01286;
713 energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5])))) * (0.964 + exp(-3.132-0.435*energy*2.0));
720 //Based on comparing MC truth information to the reconstructed energy of clusters.
721 //fNonLinearityParams[0] = 1.0;
722 //fNonLinearityParams[1] = 6.64778e-02;
723 //fNonLinearityParams[2] = 1.570;
724 //fNonLinearityParams[3] = 9.67998e-02;
725 //fNonLinearityParams[4] = 2.19381e+02;
726 //fNonLinearityParams[5] = 6.31604e+01;
727 //fNonLinearityParams[6] = 1.01286;
728 energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
734 AliDebug(2,"No correction on the energy\n");
742 //__________________________________________________
743 void AliEMCALRecoUtils::InitNonLinearityParam()
745 //Initialising Non Linearity Parameters
747 if (fNonLinearityFunction == kPi0MC) {
748 fNonLinearityParams[0] = 1.014;
749 fNonLinearityParams[1] = -0.03329;
750 fNonLinearityParams[2] = -0.3853;
751 fNonLinearityParams[3] = 0.5423;
752 fNonLinearityParams[4] = -0.4335;
755 if (fNonLinearityFunction == kPi0MCv2) {
756 fNonLinearityParams[0] = 3.11111e-02;
757 fNonLinearityParams[1] =-5.71666e-02;
758 fNonLinearityParams[2] = 5.67995e-01;
761 if (fNonLinearityFunction == kPi0MCv3) {
762 fNonLinearityParams[0] = 9.81039e-01;
763 fNonLinearityParams[1] = 1.13508e-01;
764 fNonLinearityParams[2] = 1.00173e+00;
765 fNonLinearityParams[3] = 9.67998e-02;
766 fNonLinearityParams[4] = 2.19381e+02;
767 fNonLinearityParams[5] = 6.31604e+01;
768 fNonLinearityParams[6] = 1;
771 if (fNonLinearityFunction == kPi0GammaGamma) {
772 fNonLinearityParams[0] = 1.04;
773 fNonLinearityParams[1] = -0.1445;
774 fNonLinearityParams[2] = 1.046;
777 if (fNonLinearityFunction == kPi0GammaConversion) {
778 fNonLinearityParams[0] = 0.139393;
779 fNonLinearityParams[1] = 0.0566186;
780 fNonLinearityParams[2] = 0.982133;
783 if (fNonLinearityFunction == kBeamTest) {
784 if (fNonLinearThreshold == 30) {
785 fNonLinearityParams[0] = 1.007;
786 fNonLinearityParams[1] = 0.894;
787 fNonLinearityParams[2] = 0.246;
789 if (fNonLinearThreshold == 45) {
790 fNonLinearityParams[0] = 1.003;
791 fNonLinearityParams[1] = 0.719;
792 fNonLinearityParams[2] = 0.334;
794 if (fNonLinearThreshold == 75) {
795 fNonLinearityParams[0] = 1.002;
796 fNonLinearityParams[1] = 0.797;
797 fNonLinearityParams[2] = 0.358;
801 if (fNonLinearityFunction == kBeamTestCorrected) {
802 fNonLinearityParams[0] = 0.99078;
803 fNonLinearityParams[1] = 0.161499;
804 fNonLinearityParams[2] = 0.655166;
805 fNonLinearityParams[3] = 0.134101;
806 fNonLinearityParams[4] = 163.282;
807 fNonLinearityParams[5] = 23.6904;
808 fNonLinearityParams[6] = 0.978;
811 if (fNonLinearityFunction == kBeamTestCorrectedv2) {
812 fNonLinearityParams[0] = 0.983504;
813 fNonLinearityParams[1] = 0.210106;
814 fNonLinearityParams[2] = 0.897274;
815 fNonLinearityParams[3] = 0.0829064;
816 fNonLinearityParams[4] = 152.299;
817 fNonLinearityParams[5] = 31.5028;
818 fNonLinearityParams[6] = 0.968;
821 if (fNonLinearityFunction == kSDMv5) {
822 fNonLinearityParams[0] = 1.0;
823 fNonLinearityParams[1] = 6.64778e-02;
824 fNonLinearityParams[2] = 1.570;
825 fNonLinearityParams[3] = 9.67998e-02;
826 fNonLinearityParams[4] = 2.19381e+02;
827 fNonLinearityParams[5] = 6.31604e+01;
828 fNonLinearityParams[6] = 1.01286;
831 if (fNonLinearityFunction == kPi0MCv5) {
832 fNonLinearityParams[0] = 1.0;
833 fNonLinearityParams[1] = 6.64778e-02;
834 fNonLinearityParams[2] = 1.570;
835 fNonLinearityParams[3] = 9.67998e-02;
836 fNonLinearityParams[4] = 2.19381e+02;
837 fNonLinearityParams[5] = 6.31604e+01;
838 fNonLinearityParams[6] = 1.01286;
843 //_________________________________________________________
844 Float_t AliEMCALRecoUtils::GetDepth(Float_t energy,
848 //Calculate shower depth for a given cluster energy and particle type
854 Float_t arg = energy*1000/ ecr; //Multiply energy by 1000 to transform to MeV
862 depth = x0 * (TMath::Log(arg) + 0.5);
869 depth = x0 * (TMath::Log(arg) - 0.5);
876 gGeoManager->cd("ALIC_1/XEN1_1");
877 TGeoNode *geoXEn1 = gGeoManager->GetCurrentNode();
878 TGeoNodeMatrix *geoSM = dynamic_cast<TGeoNodeMatrix *>(geoXEn1->GetDaughter(iSM));
880 TGeoVolume *geoSMVol = geoSM->GetVolume();
881 TGeoShape *geoSMShape = geoSMVol->GetShape();
882 TGeoBBox *geoBox = dynamic_cast<TGeoBBox *>(geoSMShape);
883 if (geoBox) depth = 0.5 * geoBox->GetDX()*2 ;
884 else AliFatal("Null GEANT box");
886 else AliFatal("NULL GEANT node matrix");
893 depth = x0 * (TMath::Log(arg) - 0.5);
902 depth = x0 * (TMath::Log(arg) + 0.5);
908 //____________________________________________________________________
909 void AliEMCALRecoUtils::GetMaxEnergyCell(const AliEMCALGeometry *geom,
910 AliVCaloCells* cells,
911 const AliVCluster* clu,
918 //For a given CaloCluster gets the absId of the cell
919 //with maximum energy deposit.
922 Double_t eCell = -1.;
923 Float_t fraction = 1.;
924 Float_t recalFactor = 1.;
925 Int_t cellAbsId = -1 ;
933 AliInfo("Cluster pointer null!");
934 absId=-1; iSupMod0=-1, ieta = -1; iphi = -1; shared = -1;
938 for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++) {
939 cellAbsId = clu->GetCellAbsId(iDig);
940 fraction = clu->GetCellAmplitudeFraction(iDig);
941 //printf("a Cell %d, id, %d, amp %f, fraction %f\n",iDig,cellAbsId,cells->GetCellAmplitude(cellAbsId),fraction);
942 if (fraction < 1e-4) fraction = 1.; // in case unfolding is off
943 geom->GetCellIndex(cellAbsId,iSupMod,iTower,iIphi,iIeta);
944 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);
947 } else if (iSupMod0!=iSupMod) {
949 //printf("AliEMCALRecoUtils::GetMaxEnergyCell() - SHARED CLUSTER\n");
951 if (!fCellsRecalibrated && IsRecalibrationOn()) {
952 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
954 eCell = cells->GetCellAmplitude(cellAbsId)*fraction*recalFactor;
955 //printf("b Cell %d, id, %d, amp %f, fraction %f\n",iDig,cellAbsId,eCell,fraction);
959 //printf("\t new max: cell %d, e %f, ecell %f\n",maxId, eMax,eCell);
963 //Get from the absid the supermodule, tower and eta/phi numbers
964 geom->GetCellIndex(absId,iSupMod,iTower,iIphi,iIeta);
965 //Gives SuperModule and Tower numbers
966 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,
967 iIphi, iIeta,iphi,ieta);
968 //printf("Max id %d, iSM %d, col %d, row %d\n",absId,iSupMod,ieta,iphi);
969 //printf("Max end---\n");
972 //______________________________________
973 void AliEMCALRecoUtils::InitParameters()
975 // Initialize data members with default values
977 fParticleType = kPhoton;
978 fPosAlgo = kUnchanged;
981 fNonLinearityFunction = kNoCorrection;
982 fNonLinearThreshold = 30;
984 fExoticCellFraction = 0.97;
985 fExoticCellDiffTime = 1e6;
986 fExoticCellMinAmplitude = 0.5;
988 fAODFilterMask = 128;
989 fAODHybridTracks = kFALSE;
990 fAODTPCOnlyTracks = kTRUE;
992 fCutEtaPhiSum = kTRUE;
993 fCutEtaPhiSeparate = kFALSE;
999 fClusterWindow = 100;
1004 fTrackCutsType = kLooseCut;
1007 fCutMinNClusterTPC = -1;
1008 fCutMinNClusterITS = -1;
1010 fCutMaxChi2PerClusterTPC = 1e10;
1011 fCutMaxChi2PerClusterITS = 1e10;
1013 fCutRequireTPCRefit = kFALSE;
1014 fCutRequireITSRefit = kFALSE;
1015 fCutAcceptKinkDaughters = kFALSE;
1017 fCutMaxDCAToVertexXY = 1e10;
1018 fCutMaxDCAToVertexZ = 1e10;
1019 fCutDCAToVertex2D = kFALSE;
1021 fCutRequireITSStandAlone = kFALSE; //Marcel
1022 fCutRequireITSpureSA = kFALSE; //Marcel
1024 fUseMassForTracking = kFALSE;
1026 //Misalignment matrices
1027 for (Int_t i = 0; i < 15 ; i++)
1029 fMisalTransShift[i] = 0.;
1030 fMisalRotShift[i] = 0.;
1034 for (Int_t i = 0; i < 7 ; i++) fNonLinearityParams[i] = 0.;
1036 //For kBeamTestCorrectedv2 case, but default is no correction
1037 fNonLinearityParams[0] = 0.983504;
1038 fNonLinearityParams[1] = 0.210106;
1039 fNonLinearityParams[2] = 0.897274;
1040 fNonLinearityParams[3] = 0.0829064;
1041 fNonLinearityParams[4] = 152.299;
1042 fNonLinearityParams[5] = 31.5028;
1043 fNonLinearityParams[6] = 0.968;
1045 //Cluster energy smearing
1046 fSmearClusterEnergy = kFALSE;
1047 fSmearClusterParam[0] = 0.07; // * sqrt E term
1048 fSmearClusterParam[1] = 0.00; // * E term
1049 fSmearClusterParam[2] = 0.00; // constant
1052 //_____________________________________________________
1053 void AliEMCALRecoUtils::InitEMCALRecalibrationFactors()
1055 //Init EMCAL recalibration factors
1056 AliDebug(2,"AliCalorimeterUtils::InitEMCALRecalibrationFactors()");
1057 //In order to avoid rewriting the same histograms
1058 Bool_t oldStatus = TH1::AddDirectoryStatus();
1059 TH1::AddDirectory(kFALSE);
1061 fEMCALRecalibrationFactors = new TObjArray(12);
1062 for (int i = 0; i < 12; i++)
1063 fEMCALRecalibrationFactors->Add(new TH2F(Form("EMCALRecalFactors_SM%d",i),
1064 Form("EMCALRecalFactors_SM%d",i), 48, 0, 48, 24, 0, 24));
1065 //Init the histograms with 1
1066 for (Int_t sm = 0; sm < 12; sm++)
1068 for (Int_t i = 0; i < 48; i++)
1070 for (Int_t j = 0; j < 24; j++)
1072 SetEMCALChannelRecalibrationFactor(sm,i,j,1.);
1077 fEMCALRecalibrationFactors->SetOwner(kTRUE);
1078 fEMCALRecalibrationFactors->Compress();
1080 //In order to avoid rewriting the same histograms
1081 TH1::AddDirectory(oldStatus);
1084 //_________________________________________________________
1085 void AliEMCALRecoUtils::InitEMCALTimeRecalibrationFactors()
1087 //Init EMCAL recalibration factors
1088 AliDebug(2,"AliCalorimeterUtils::InitEMCALRecalibrationFactors()");
1089 //In order to avoid rewriting the same histograms
1090 Bool_t oldStatus = TH1::AddDirectoryStatus();
1091 TH1::AddDirectory(kFALSE);
1093 fEMCALTimeRecalibrationFactors = new TObjArray(4);
1094 for (int i = 0; i < 4; i++)
1095 fEMCALTimeRecalibrationFactors->Add(new TH1F(Form("hAllTimeAvBC%d",i),
1096 Form("hAllTimeAvBC%d",i),
1097 48*24*12,0.,48*24*12) );
1098 //Init the histograms with 1
1099 for (Int_t bc = 0; bc < 4; bc++)
1101 for (Int_t i = 0; i < 48*24*12; i++)
1102 SetEMCALChannelTimeRecalibrationFactor(bc,i,0.);
1105 fEMCALTimeRecalibrationFactors->SetOwner(kTRUE);
1106 fEMCALTimeRecalibrationFactors->Compress();
1108 //In order to avoid rewriting the same histograms
1109 TH1::AddDirectory(oldStatus);
1112 //____________________________________________________
1113 void AliEMCALRecoUtils::InitEMCALBadChannelStatusMap()
1115 //Init EMCAL bad channels map
1116 AliDebug(2,"AliEMCALRecoUtils::InitEMCALBadChannelStatusMap()");
1117 //In order to avoid rewriting the same histograms
1118 Bool_t oldStatus = TH1::AddDirectoryStatus();
1119 TH1::AddDirectory(kFALSE);
1121 fEMCALBadChannelMap = new TObjArray(12);
1122 //TH2F * hTemp = new TH2I("EMCALBadChannelMap","EMCAL SuperModule bad channel map", 48, 0, 48, 24, 0, 24);
1123 for (int i = 0; i < 12; i++)
1125 fEMCALBadChannelMap->Add(new TH2I(Form("EMCALBadChannelMap_Mod%d",i),Form("EMCALBadChannelMap_Mod%d",i), 48, 0, 48, 24, 0, 24));
1128 fEMCALBadChannelMap->SetOwner(kTRUE);
1129 fEMCALBadChannelMap->Compress();
1131 //In order to avoid rewriting the same histograms
1132 TH1::AddDirectory(oldStatus);
1135 //____________________________________________________________________________
1136 void AliEMCALRecoUtils::RecalibrateClusterEnergy(const AliEMCALGeometry* geom,
1137 AliVCluster * cluster,
1138 AliVCaloCells * cells,
1141 // Recalibrate the cluster energy and Time, considering the recalibration map
1142 // and the energy of the cells and time that compose the cluster.
1143 // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1146 AliInfo("Cluster pointer null!");
1150 //Get the cluster number of cells and list of absId, check what kind of cluster do we have.
1151 UShort_t * index = cluster->GetCellsAbsId() ;
1152 Double_t * fraction = cluster->GetCellsAmplitudeFraction() ;
1153 Int_t ncells = cluster->GetNCells();
1155 //Initialize some used variables
1158 Int_t icol =-1, irow =-1, imod=1;
1159 Float_t factor = 1, frac = 0;
1160 Int_t absIdMax = -1;
1163 //Loop on the cells, get the cell amplitude and recalibration factor, multiply and and to the new energy
1164 for (Int_t icell = 0; icell < ncells; icell++)
1166 absId = index[icell];
1167 frac = fraction[icell];
1168 if (frac < 1e-5) frac = 1; //in case of EMCAL, this is set as 0 since unfolding is off
1170 if (!fCellsRecalibrated && IsRecalibrationOn()) {
1172 Int_t iTower = -1, iIphi = -1, iIeta = -1;
1173 geom->GetCellIndex(absId,imod,iTower,iIphi,iIeta);
1174 if (fEMCALRecalibrationFactors->GetEntries() <= imod)
1176 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,irow,icol);
1177 factor = GetEMCALChannelRecalibrationFactor(imod,icol,irow);
1179 AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - recalibrate cell: module %d, col %d, row %d, cell fraction %f,recalibration factor %f, cell energy %f\n",
1180 imod,icol,irow,frac,factor,cells->GetCellAmplitude(absId)));
1184 energy += cells->GetCellAmplitude(absId)*factor*frac;
1186 if (emax < cells->GetCellAmplitude(absId)*factor*frac) {
1187 emax = cells->GetCellAmplitude(absId)*factor*frac;
1192 AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - Energy before %f, after %f \n",cluster->E(),energy));
1194 cluster->SetE(energy);
1196 // Recalculate time of cluster
1197 Double_t timeorg = cluster->GetTOF();
1199 Double_t time = cells->GetCellTime(absIdMax);
1200 if (!fCellsRecalibrated && IsTimeRecalibrationOn())
1201 RecalibrateCellTime(absIdMax,bc,time);
1203 cluster->SetTOF(time);
1205 AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - Time before %f, after %f \n",timeorg,cluster->GetTOF()));
1208 //_____________________________________________________________
1209 void AliEMCALRecoUtils::RecalibrateCells(AliVCaloCells * cells,
1212 // Recalibrate the cells time and energy, considering the recalibration map and the energy
1213 // of the cells that compose the cluster.
1214 // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1216 if (!IsRecalibrationOn() && !IsTimeRecalibrationOn() && !IsBadChannelsRemovalSwitchedOn())
1220 AliInfo("Cells pointer null!");
1225 Bool_t accept = kFALSE;
1228 Double_t ecellin = 0;
1229 Double_t tcellin = 0;
1233 Int_t nEMcell = cells->GetNumberOfCells() ;
1234 for (Int_t iCell = 0; iCell < nEMcell; iCell++)
1236 cells->GetCell( iCell, absId, ecellin, tcellin, mclabel, efrac );
1238 accept = AcceptCalibrateCell(absId, bc, ecell ,tcell ,cells);
1245 cells->SetCell(iCell,absId,ecell, tcell, mclabel, efrac);
1248 fCellsRecalibrated = kTRUE;
1251 //_______________________________________________________________________________________________________
1252 void AliEMCALRecoUtils::RecalibrateCellTime(Int_t absId, Int_t bc, Double_t & celltime) const
1254 // Recalibrate time of cell with absID considering the recalibration map
1255 // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1257 if (!fCellsRecalibrated && IsTimeRecalibrationOn() && bc >= 0) {
1258 celltime -= GetEMCALChannelTimeRecalibrationFactor(bc%4,absId)*1.e-9; ;
1262 //______________________________________________________________________________
1263 void AliEMCALRecoUtils::RecalculateClusterPosition(const AliEMCALGeometry *geom,
1264 AliVCaloCells* cells,
1267 //For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1270 AliInfo("Cluster pointer null!");
1274 if (fPosAlgo==kPosTowerGlobal) RecalculateClusterPositionFromTowerGlobal( geom, cells, clu);
1275 else if (fPosAlgo==kPosTowerIndex) RecalculateClusterPositionFromTowerIndex ( geom, cells, clu);
1276 else AliDebug(2,"Algorithm to recalculate position not selected, do nothing.");
1279 //_____________________________________________________________________________________________
1280 void AliEMCALRecoUtils::RecalculateClusterPositionFromTowerGlobal(const AliEMCALGeometry *geom,
1281 AliVCaloCells* cells,
1284 // For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1285 // The algorithm is a copy of what is done in AliEMCALRecPoint
1287 Double_t eCell = 0.;
1288 Float_t fraction = 1.;
1289 Float_t recalFactor = 1.;
1292 Int_t iTower = -1, iIphi = -1, iIeta = -1;
1293 Int_t iSupModMax = -1, iSM=-1, iphi = -1, ieta = -1;
1294 Float_t weight = 0., totalWeight=0.;
1295 Float_t newPos[3] = {0,0,0};
1296 Double_t pLocal[3], pGlobal[3];
1297 Bool_t shared = kFALSE;
1299 Float_t clEnergy = clu->E(); //Energy already recalibrated previously
1302 GetMaxEnergyCell(geom, cells, clu, absId, iSupModMax, ieta, iphi,shared);
1303 Double_t depth = GetDepth(clEnergy,fParticleType,iSupModMax) ;
1305 //printf("** Cluster energy %f, ncells %d, depth %f\n",clEnergy,clu->GetNCells(),depth);
1307 for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++)
1309 absId = clu->GetCellAbsId(iDig);
1310 fraction = clu->GetCellAmplitudeFraction(iDig);
1311 if (fraction < 1e-4) fraction = 1.; // in case unfolding is off
1313 if (!fCellsRecalibrated) {
1314 geom->GetCellIndex(absId,iSM,iTower,iIphi,iIeta);
1315 geom->GetCellPhiEtaIndexInSModule(iSM,iTower,iIphi, iIeta,iphi,ieta);
1316 if (IsRecalibrationOn()) {
1317 recalFactor = GetEMCALChannelRecalibrationFactor(iSM,ieta,iphi);
1321 eCell = cells->GetCellAmplitude(absId)*fraction*recalFactor;
1323 weight = GetCellWeight(eCell,clEnergy);
1324 totalWeight += weight;
1326 geom->RelPosCellInSModule(absId,depth,pLocal[0],pLocal[1],pLocal[2]);
1327 //printf("pLocal (%f,%f,%f), SM %d, absId %d\n",pLocal[0],pLocal[1],pLocal[2],iSupModMax,absId);
1328 geom->GetGlobal(pLocal,pGlobal,iSupModMax);
1329 //printf("pLocal (%f,%f,%f)\n",pGlobal[0],pGlobal[1],pGlobal[2]);
1331 for (int i=0; i<3; i++ ) newPos[i] += (weight*pGlobal[i]);
1334 if (totalWeight>0) {
1335 for (int i=0; i<3; i++ ) newPos[i] /= totalWeight;
1338 //Float_t pos[]={0,0,0};
1339 //clu->GetPosition(pos);
1340 //printf("OldPos : %2.3f,%2.3f,%2.3f\n",pos[0],pos[1],pos[2]);
1341 //printf("NewPos : %2.3f,%2.3f,%2.3f\n",newPos[0],newPos[1],newPos[2]);
1343 if (iSupModMax > 1) { //sector 1
1344 newPos[0] +=fMisalTransShift[3];//-=3.093;
1345 newPos[1] +=fMisalTransShift[4];//+=6.82;
1346 newPos[2] +=fMisalTransShift[5];//+=1.635;
1347 //printf(" + : %2.3f,%2.3f,%2.3f\n",fMisalTransShift[3],fMisalTransShift[4],fMisalTransShift[5]);
1349 newPos[0] +=fMisalTransShift[0];//+=1.134;
1350 newPos[1] +=fMisalTransShift[1];//+=8.2;
1351 newPos[2] +=fMisalTransShift[2];//+=1.197;
1352 //printf(" + : %2.3f,%2.3f,%2.3f\n",fMisalTransShift[0],fMisalTransShift[1],fMisalTransShift[2]);
1354 //printf("NewPos : %2.3f,%2.3f,%2.3f\n",newPos[0],newPos[1],newPos[2]);
1356 clu->SetPosition(newPos);
1359 //____________________________________________________________________________________________
1360 void AliEMCALRecoUtils::RecalculateClusterPositionFromTowerIndex(const AliEMCALGeometry *geom,
1361 AliVCaloCells* cells,
1364 // For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1365 // The algorithm works with the tower indeces, averages the indeces and from them it calculates the global position
1367 Double_t eCell = 1.;
1368 Float_t fraction = 1.;
1369 Float_t recalFactor = 1.;
1373 Int_t iIphi = -1, iIeta = -1;
1374 Int_t iSupMod = -1, iSupModMax = -1;
1375 Int_t iphi = -1, ieta =-1;
1376 Bool_t shared = kFALSE;
1378 Float_t clEnergy = clu->E(); //Energy already recalibrated previously.
1382 GetMaxEnergyCell(geom, cells, clu, absId, iSupModMax, ieta, iphi,shared);
1383 Float_t depth = GetDepth(clEnergy,fParticleType,iSupMod) ;
1385 Float_t weight = 0., weightedCol = 0., weightedRow = 0., totalWeight=0.;
1386 Bool_t areInSameSM = kTRUE; //exclude clusters with cells in different SMs for now
1387 Int_t startingSM = -1;
1389 for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++)
1391 absId = clu->GetCellAbsId(iDig);
1392 fraction = clu->GetCellAmplitudeFraction(iDig);
1393 if (fraction < 1e-4) fraction = 1.; // in case unfolding is off
1395 if (iDig==0) startingSM = iSupMod;
1396 else if (iSupMod != startingSM) areInSameSM = kFALSE;
1398 eCell = cells->GetCellAmplitude(absId);
1400 geom->GetCellIndex(absId,iSupMod,iTower,iIphi,iIeta);
1401 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);
1403 if (!fCellsRecalibrated)
1405 if (IsRecalibrationOn()) {
1406 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1410 eCell = cells->GetCellAmplitude(absId)*fraction*recalFactor;
1412 weight = GetCellWeight(eCell,clEnergy);
1413 if (weight < 0) weight = 0;
1414 totalWeight += weight;
1415 weightedCol += ieta*weight;
1416 weightedRow += iphi*weight;
1418 //printf("Max cell? cell %d, amplitude org %f, fraction %f, recalibration %f, amplitude new %f \n",cellAbsId, cells->GetCellAmplitude(cellAbsId), fraction, recalFactor, eCell) ;
1421 Float_t xyzNew[]={0.,0.,0.};
1422 if (areInSameSM == kTRUE) {
1423 //printf("In Same SM\n");
1424 weightedCol = weightedCol/totalWeight;
1425 weightedRow = weightedRow/totalWeight;
1426 geom->RecalculateTowerPosition(weightedRow, weightedCol, iSupModMax, depth, fMisalTransShift, fMisalRotShift, xyzNew);
1430 //printf("In Different SM\n");
1431 geom->RecalculateTowerPosition(iphi, ieta, iSupModMax, depth, fMisalTransShift, fMisalRotShift, xyzNew);
1434 clu->SetPosition(xyzNew);
1437 //___________________________________________________________________________________________
1438 void AliEMCALRecoUtils::RecalculateClusterDistanceToBadChannel(const AliEMCALGeometry * geom,
1439 AliVCaloCells* cells,
1440 AliVCluster * cluster)
1442 //re-evaluate distance to bad channel with updated bad map
1444 if (!fRecalDistToBadChannels) return;
1448 AliInfo("Cluster pointer null!");
1452 //Get channels map of the supermodule where the cluster is.
1453 Int_t absIdMax = -1, iSupMod =-1, icolM = -1, irowM = -1;
1454 Bool_t shared = kFALSE;
1455 GetMaxEnergyCell(geom, cells, cluster, absIdMax, iSupMod, icolM, irowM, shared);
1456 TH2D* hMap = (TH2D*)fEMCALBadChannelMap->At(iSupMod);
1459 Float_t minDist = 10000.;
1462 //Loop on tower status map
1463 for (Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++)
1465 for (Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++)
1467 //Check if tower is bad.
1468 if (hMap->GetBinContent(icol,irow)==0) continue;
1469 //printf("AliEMCALRecoUtils::RecalculateDistanceToBadChannels() - \n \t Bad channel in SM %d, col %d, row %d, \n \t Cluster max in col %d, row %d\n",
1470 // iSupMod,icol, irow, icolM,irowM);
1472 dRrow=TMath::Abs(irowM-irow);
1473 dRcol=TMath::Abs(icolM-icol);
1474 dist=TMath::Sqrt(dRrow*dRrow+dRcol*dRcol);
1477 //printf("MIN DISTANCE TO BAD %2.2f\n",dist);
1483 //In case the cluster is shared by 2 SuperModules, need to check the map of the second Super Module
1487 Int_t iSupMod2 = -1;
1489 //The only possible combinations are (0,1), (2,3) ... (8,9)
1490 if (iSupMod%2) iSupMod2 = iSupMod-1;
1491 else iSupMod2 = iSupMod+1;
1492 hMap2 = (TH2D*)fEMCALBadChannelMap->At(iSupMod2);
1494 //Loop on tower status map of second super module
1495 for (Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++)
1497 for (Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++)
1499 //Check if tower is bad.
1500 if (hMap2->GetBinContent(icol,irow)==0)
1502 //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",
1503 // iSupMod2,icol, irow,iSupMod,icolM,irowM);
1504 dRrow=TMath::Abs(irow-irowM);
1507 dRcol=TMath::Abs(icol-(AliEMCALGeoParams::fgkEMCALCols+icolM));
1509 dRcol=TMath::Abs(AliEMCALGeoParams::fgkEMCALCols+icol-icolM);
1512 dist=TMath::Sqrt(dRrow*dRrow+dRcol*dRcol);
1513 if (dist < minDist) minDist = dist;
1516 }// shared cluster in 2 SuperModules
1518 AliDebug(2,Form("Max cluster cell (SM,col,row)=(%d %d %d) - Distance to Bad Channel %2.2f",iSupMod, icolM, irowM, minDist));
1519 cluster->SetDistanceToBadChannel(minDist);
1522 //__________________________________________________________________
1523 void AliEMCALRecoUtils::RecalculateClusterPID(AliVCluster * cluster)
1525 //re-evaluate identification parameters with bayesian
1528 AliInfo("Cluster pointer null!");
1532 if (cluster->GetM02() != 0)
1533 fPIDUtils->ComputePID(cluster->E(),cluster->GetM02());
1535 Float_t pidlist[AliPID::kSPECIESCN+1];
1536 for (Int_t i = 0; i < AliPID::kSPECIESCN+1; i++) pidlist[i] = fPIDUtils->GetPIDFinal(i);
1538 cluster->SetPID(pidlist);
1541 //___________________________________________________________________________________________________________________
1542 void AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(const AliEMCALGeometry * geom,
1543 AliVCaloCells* cells,
1544 AliVCluster * cluster,
1545 Float_t & l0, Float_t & l1,
1546 Float_t & disp, Float_t & dEta, Float_t & dPhi,
1547 Float_t & sEta, Float_t & sPhi, Float_t & sEtaPhi)
1549 // Calculates new center of gravity in the local EMCAL-module coordinates
1550 // and tranfers into global ALICE coordinates
1551 // Calculates Dispersion and main axis
1554 AliInfo("Cluster pointer null!");
1558 Double_t eCell = 0.;
1559 Float_t fraction = 1.;
1560 Float_t recalFactor = 1.;
1568 Double_t etai = -1.;
1569 Double_t phii = -1.;
1574 Double_t etaMean = 0.;
1575 Double_t phiMean = 0.;
1577 //Loop on cells, calculate the cluster energy, in case a cut on cell energy is added
1578 // and to check if the cluster is between 2 SM in eta
1580 Bool_t shared = kFALSE;
1583 for (Int_t iDigit=0; iDigit < cluster->GetNCells(); iDigit++)
1585 //Get from the absid the supermodule, tower and eta/phi numbers
1586 geom->GetCellIndex(cluster->GetCellAbsId(iDigit),iSupMod,iTower,iIphi,iIeta);
1587 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi,iIeta, iphi,ieta);
1589 //Check if there are cells of different SM
1590 if (iDigit == 0 ) iSM0 = iSupMod;
1591 else if (iSupMod!= iSM0) shared = kTRUE;
1593 //Get the cell energy, if recalibration is on, apply factors
1594 fraction = cluster->GetCellAmplitudeFraction(iDigit);
1595 if (fraction < 1e-4) fraction = 1.; // in case unfolding is off
1597 if (IsRecalibrationOn()) {
1598 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1601 eCell = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1608 for (Int_t iDigit=0; iDigit < cluster->GetNCells(); iDigit++)
1610 //Get from the absid the supermodule, tower and eta/phi numbers
1611 geom->GetCellIndex(cluster->GetCellAbsId(iDigit),iSupMod,iTower,iIphi,iIeta);
1612 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi,iIeta, iphi,ieta);
1614 //Get the cell energy, if recalibration is on, apply factors
1615 fraction = cluster->GetCellAmplitudeFraction(iDigit);
1616 if (fraction < 1e-4) fraction = 1.; // in case unfolding is off
1618 if (!fCellsRecalibrated) {
1619 if (IsRecalibrationOn()) {
1620 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1624 eCell = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1626 // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
1627 // C Side impair SM, nSupMod%2=1; A side pair SM, nSupMod%2=0
1628 if (shared && iSupMod%2) ieta+=AliEMCALGeoParams::fgkEMCALCols;
1630 if (cluster->E() > 0 && eCell > 0) {
1631 w = GetCellWeight(eCell,cluster->E());
1633 etai=(Double_t)ieta;
1634 phii=(Double_t)iphi;
1640 sEta += w * etai * etai ;
1641 etaMean += w * etai ;
1642 sPhi += w * phii * phii ;
1643 phiMean += w * phii ;
1644 sEtaPhi += w * etai * phii ;
1647 AliError(Form("Wrong energy %f and/or amplitude %f\n", eCell, cluster->E()));
1650 //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()) {
1668 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1670 eCell = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1672 // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
1673 // C Side impair SM, nSupMod%2=1; A side pair SM, nSupMod%2=0
1674 if (shared && iSupMod%2) ieta+=AliEMCALGeoParams::fgkEMCALCols;
1676 if (cluster->E() > 0 && eCell > 0) {
1677 w = GetCellWeight(eCell,cluster->E());
1679 etai=(Double_t)ieta;
1680 phii=(Double_t)iphi;
1682 disp += w *((etai-etaMean)*(etai-etaMean)+(phii-phiMean)*(phii-phiMean));
1683 dEta += w * (etai-etaMean)*(etai-etaMean) ;
1684 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) {
1699 sEta -= etaMean * etaMean ;
1700 sPhi -= phiMean * phiMean ;
1701 sEtaPhi -= etaMean * phiMean ;
1703 l0 = (0.5 * (sEta + sPhi) + TMath::Sqrt( 0.25 * (sEta - sPhi) * (sEta - sPhi) + sEtaPhi * sEtaPhi ));
1704 l1 = (0.5 * (sEta + sPhi) - TMath::Sqrt( 0.25 * (sEta - sPhi) * (sEta - sPhi) + sEtaPhi * sEtaPhi ));
1708 dEta = 0. ; dPhi = 0. ; disp = 0. ;
1709 sEta = 0. ; sPhi = 0. ; sEtaPhi = 0. ;
1713 //____________________________________________________________________________________________
1714 void AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(const AliEMCALGeometry * geom,
1715 AliVCaloCells* cells,
1716 AliVCluster * cluster)
1718 // Calculates new center of gravity in the local EMCAL-module coordinates
1719 // and tranfers into global ALICE coordinates
1720 // Calculates Dispersion and main axis and puts them into the cluster
1722 Float_t l0 = 0., l1 = 0.;
1723 Float_t disp = 0., dEta = 0., dPhi = 0.;
1724 Float_t sEta = 0., sPhi = 0., sEtaPhi = 0.;
1726 AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(geom,cells,cluster,l0,l1,disp,
1727 dEta, dPhi, sEta, sPhi, sEtaPhi);
1729 cluster->SetM02(l0);
1730 cluster->SetM20(l1);
1731 if (disp > 0. ) cluster->SetDispersion(TMath::Sqrt(disp)) ;
1735 //____________________________________________________________________________
1736 void AliEMCALRecoUtils::FindMatches(AliVEvent *event,
1737 TObjArray * clusterArr,
1738 const AliEMCALGeometry *geom)
1740 //This function should be called before the cluster loop
1741 //Before call this function, please recalculate the cluster positions
1742 //Given the input event, loop over all the tracks, select the closest cluster as matched with fCutR
1743 //Store matched cluster indexes and residuals
1745 fMatchedTrackIndex ->Reset();
1746 fMatchedClusterIndex->Reset();
1747 fResidualPhi->Reset();
1748 fResidualEta->Reset();
1750 fMatchedTrackIndex ->Set(1000);
1751 fMatchedClusterIndex->Set(1000);
1752 fResidualPhi->Set(1000);
1753 fResidualEta->Set(1000);
1755 AliESDEvent* esdevent = dynamic_cast<AliESDEvent*> (event);
1756 AliAODEvent* aodevent = dynamic_cast<AliAODEvent*> (event);
1758 // init the magnetic field if not already on
1759 if (!TGeoGlobalMagField::Instance()->GetField()) {
1760 if (!event->InitMagneticField())
1762 AliInfo("Mag Field not initialized, null esd/aod evetn pointers");
1767 UInt_t mask1 = esdevent->GetESDRun()->GetDetectorsInDAQ();
1768 UInt_t mask2 = esdevent->GetESDRun()->GetDetectorsInReco();
1769 Bool_t desc1 = (mask1 >> 3) & 0x1;
1770 Bool_t desc2 = (mask2 >> 3) & 0x1;
1771 if (desc1==0 || desc2==0) {
1772 // AliError(Form("TPC not in DAQ/RECO: %u (%u)/%u (%u)",
1773 // mask1, esdevent->GetESDRun()->GetDetectorsInReco(),
1774 // mask2, esdevent->GetESDRun()->GetDetectorsInDAQ()));
1779 TObjArray *clusterArray = 0x0;
1781 clusterArray = new TObjArray(event->GetNumberOfCaloClusters());
1782 for (Int_t icl=0; icl<event->GetNumberOfCaloClusters(); icl++)
1784 AliVCluster *cluster = (AliVCluster*) event->GetCaloCluster(icl);
1785 if (geom && !IsGoodCluster(cluster,geom,(AliVCaloCells*)event->GetEMCALCells()))
1787 clusterArray->AddAt(cluster,icl);
1793 for (Int_t i=0; i<21;i++) cv[i]=0;
1794 for (Int_t itr=0; itr<event->GetNumberOfTracks(); itr++)
1796 AliExternalTrackParam *trackParam = 0;
1798 //If the input event is ESD, the starting point for extrapolation is TPCOut, if available, or TPCInner
1799 AliESDtrack *esdTrack = 0;
1800 AliAODTrack *aodTrack = 0;
1802 esdTrack = esdevent->GetTrack(itr);
1803 if (!esdTrack) continue;
1804 if (!IsAccepted(esdTrack)) continue;
1805 if (esdTrack->Pt()<fCutMinTrackPt) continue;
1806 Double_t phi = esdTrack->Phi()*TMath::RadToDeg();
1807 if (TMath::Abs(esdTrack->Eta())>0.9 || phi <= 10 || phi >= 250 ) continue;
1809 trackParam = const_cast<AliExternalTrackParam*>(esdTrack->GetInnerParam()); // if TPC Available
1811 trackParam = new AliExternalTrackParam(*esdTrack); // If ITS Track Standing alone
1814 //If the input event is AOD, the starting point for extrapolation is at vertex
1815 //AOD tracks are selected according to its filterbit.
1816 else if (aodevent) {
1817 aodTrack = dynamic_cast<AliAODTrack*>(aodevent->GetTrack(itr));
1818 if(!aodTrack) AliFatal("Not a standard AOD");
1819 if (!aodTrack) continue;
1821 if (fAODTPCOnlyTracks) { // Match with TPC only tracks, default from May 2013, before filter bit 32
1822 //printf("Match with TPC only tracks, accept? %d, test bit 128 <%d> \n", aodTrack->IsTPCOnly(), aodTrack->TestFilterMask(128));
1823 if (!aodTrack->IsTPCConstrained()) continue ;
1824 } else if (fAODHybridTracks) { // Match with hybrid tracks
1825 //printf("Match with Hybrid tracks, accept? %d \n", aodTrack->IsHybridGlobalConstrainedGlobal());
1826 if (!aodTrack->IsHybridGlobalConstrainedGlobal()) continue ;
1827 } else { // Match with tracks on a mask
1828 //printf("Match with tracks having filter bit mask %d, accept? %d \n",fAODFilterMask,aodTrack->TestFilterMask(fAODFilterMask));
1829 if (!aodTrack->TestFilterMask(fAODFilterMask) ) continue; //Select AOD tracks
1832 if (aodTrack->Pt()<fCutMinTrackPt) continue;
1834 Double_t phi = aodTrack->Phi()*TMath::RadToDeg();
1835 if (TMath::Abs(aodTrack->Eta())>0.9 || phi <= 10 || phi >= 250 )
1837 Double_t pos[3],mom[3];
1838 aodTrack->GetXYZ(pos);
1839 aodTrack->GetPxPyPz(mom);
1840 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()));
1842 trackParam= new AliExternalTrackParam(pos,mom,cv,aodTrack->Charge());
1845 //Return if the input data is not "AOD" or "ESD"
1847 printf("Wrong input data type! Should be \"AOD\" or \"ESD\"\n");
1849 clusterArray->Clear();
1850 delete clusterArray;
1855 if (!trackParam) continue;
1857 //Extrapolate the track to EMCal surface
1858 AliExternalTrackParam emcalParam(*trackParam);
1859 Float_t eta, phi, pt;
1860 if (!ExtrapolateTrackToEMCalSurface(&emcalParam, fEMCalSurfaceDistance, fMass, fStepSurface, eta, phi, pt)) {
1861 if (aodevent && trackParam) delete trackParam;
1862 if (fITSTrackSA && trackParam) delete trackParam;
1866 if (TMath::Abs(eta)>0.75 || (phi) < 70*TMath::DegToRad() || (phi) > 190*TMath::DegToRad()) {
1867 if (aodevent && trackParam) delete trackParam;
1868 if (fITSTrackSA && trackParam) delete trackParam;
1872 //Find matched clusters
1874 Float_t dEta = -999, dPhi = -999;
1876 index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArray, dEta, dPhi);
1878 index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArr, dEta, dPhi);
1882 fMatchedTrackIndex ->AddAt(itr,matched);
1883 fMatchedClusterIndex ->AddAt(index,matched);
1884 fResidualEta ->AddAt(dEta,matched);
1885 fResidualPhi ->AddAt(dPhi,matched);
1888 if (aodevent && trackParam) delete trackParam;
1889 if (fITSTrackSA && trackParam) delete trackParam;
1893 clusterArray->Clear();
1894 delete clusterArray;
1897 AliDebug(2,Form("Number of matched pairs = %d !\n",matched));
1899 fMatchedTrackIndex ->Set(matched);
1900 fMatchedClusterIndex ->Set(matched);
1901 fResidualPhi ->Set(matched);
1902 fResidualEta ->Set(matched);
1905 //________________________________________________________________________________
1906 Int_t AliEMCALRecoUtils::FindMatchedClusterInEvent(const AliESDtrack *track,
1907 const AliVEvent *event,
1908 const AliEMCALGeometry *geom,
1909 Float_t &dEta, Float_t &dPhi)
1912 // This function returns the index of matched cluster to input track
1913 // Returns -1 if no match is found
1915 Double_t phiV = track->Phi()*TMath::RadToDeg();
1916 if (TMath::Abs(track->Eta())>0.9 || phiV <= 10 || phiV >= 250 ) return index;
1917 AliExternalTrackParam *trackParam = 0;
1919 trackParam = const_cast<AliExternalTrackParam*>(track->GetInnerParam()); // If TPC
1921 trackParam = new AliExternalTrackParam(*track);
1923 if (!trackParam) return index;
1924 AliExternalTrackParam emcalParam(*trackParam);
1925 Float_t eta, phi, pt;
1927 if (!ExtrapolateTrackToEMCalSurface(&emcalParam, fEMCalSurfaceDistance, fMass, fStepSurface, eta, phi, pt)) {
1928 if (fITSTrackSA) delete trackParam;
1931 if (TMath::Abs(eta)>0.75 || (phi) < 70*TMath::DegToRad() || (phi) > 190*TMath::DegToRad()) {
1932 if (fITSTrackSA) delete trackParam;
1936 TObjArray *clusterArr = new TObjArray(event->GetNumberOfCaloClusters());
1938 for (Int_t icl=0; icl<event->GetNumberOfCaloClusters(); icl++)
1940 AliVCluster *cluster = (AliVCluster*) event->GetCaloCluster(icl);
1941 if (geom && !IsGoodCluster(cluster,geom,(AliVCaloCells*)event->GetEMCALCells())) continue;
1942 clusterArr->AddAt(cluster,icl);
1945 index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArr, dEta, dPhi);
1946 clusterArr->Clear();
1948 if (fITSTrackSA) delete trackParam;
1953 //_______________________________________________________________________________________________
1954 Int_t AliEMCALRecoUtils::FindMatchedClusterInClusterArr(const AliExternalTrackParam *emcalParam,
1955 AliExternalTrackParam *trkParam,
1956 const TObjArray * clusterArr,
1957 Float_t &dEta, Float_t &dPhi)
1959 // Find matched cluster in array
1961 dEta=-999, dPhi=-999;
1962 Float_t dRMax = fCutR, dEtaMax=fCutEta, dPhiMax=fCutPhi;
1964 Float_t tmpEta=-999, tmpPhi=-999;
1966 Double_t exPos[3] = {0.,0.,0.};
1967 if (!emcalParam->GetXYZ(exPos)) return index;
1969 Float_t clsPos[3] = {0.,0.,0.};
1970 for (Int_t icl=0; icl<clusterArr->GetEntriesFast(); icl++)
1972 AliVCluster *cluster = dynamic_cast<AliVCluster*> (clusterArr->At(icl)) ;
1973 if (!cluster || !cluster->IsEMCAL()) continue;
1974 cluster->GetPosition(clsPos);
1975 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));
1976 if (dR > fClusterWindow) continue;
1978 AliExternalTrackParam trkPamTmp (*trkParam);//Retrieve the starting point every time before the extrapolation
1979 if (!ExtrapolateTrackToCluster(&trkPamTmp, cluster, fMass, fStepCluster, tmpEta, tmpPhi)) continue;
1980 if (fCutEtaPhiSum) {
1981 Float_t tmpR=TMath::Sqrt(tmpEta*tmpEta + tmpPhi*tmpPhi);
1988 } else if (fCutEtaPhiSeparate) {
1989 if (TMath::Abs(tmpEta)<TMath::Abs(dEtaMax) && TMath::Abs(tmpPhi)<TMath::Abs(dPhiMax)) {
1995 printf("Error: please specify your cut criteria\n");
1996 printf("To cut on sqrt(dEta^2+dPhi^2), use: SwitchOnCutEtaPhiSum()\n");
1997 printf("To cut on dEta and dPhi separately, use: SwitchOnCutEtaPhiSeparate()\n");
2008 //------------------------------------------------------------------------------------
2009 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToEMCalSurface(AliVTrack *track,
2010 Double_t emcalR, Double_t mass,
2011 Double_t step, Double_t minpt)
2013 // Extrapolate track to EMCAL surface
2015 track->SetTrackPhiEtaPtOnEMCal(-999, -999, -999);
2017 if (track->Pt()<minpt)
2020 Double_t phi = track->Phi()*TMath::RadToDeg();
2021 if (TMath::Abs(track->Eta())>0.9 || phi <= 10 || phi >= 250)
2024 AliESDtrack *esdt = dynamic_cast<AliESDtrack*>(track);
2025 AliAODTrack *aodt = 0;
2027 aodt = dynamic_cast<AliAODTrack*>(track);
2032 // Select the mass hypothesis
2035 Bool_t onlyTPC = kFALSE;
2036 if ( mass == -99 ) onlyTPC=kTRUE;
2040 if ( fUseMassForTracking ) mass = esdt->GetMassForTracking();
2041 else mass = esdt->GetMass(onlyTPC);
2045 if ( fUseMassForTracking ) mass = aodt->GetMassForTracking();
2046 else mass = aodt->M();
2050 AliExternalTrackParam *trackParam = 0;
2052 const AliExternalTrackParam *in = esdt->GetInnerParam();
2055 trackParam = new AliExternalTrackParam(*in);
2057 Double_t xyz[3] = {0}, pxpypz[3] = {0}, cv[21] = {0};
2058 aodt->PxPyPz(pxpypz);
2060 aodt->GetCovarianceXYZPxPyPz(cv);
2061 trackParam = new AliExternalTrackParam(xyz,pxpypz,cv,aodt->Charge());
2066 Float_t etaout=-999, phiout=-999, ptout=-999;
2067 Bool_t ret = ExtrapolateTrackToEMCalSurface(trackParam,
2077 if (TMath::Abs(etaout)>0.75 || (phiout<70*TMath::DegToRad()) || (phiout>190*TMath::DegToRad()))
2079 track->SetTrackPhiEtaPtOnEMCal(phiout, etaout, ptout);
2084 //------------------------------------------------------------------------------------
2085 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToEMCalSurface(AliExternalTrackParam *trkParam,
2093 //Extrapolate track to EMCAL surface
2095 eta = -999, phi = -999, pt = -999;
2096 if (!trkParam) return kFALSE;
2097 if (!AliTrackerBase::PropagateTrackToBxByBz(trkParam, emcalR, mass, step, kTRUE, 0.8, -1)) return kFALSE;
2098 Double_t trkPos[3] = {0.,0.,0.};
2099 if (!trkParam->GetXYZ(trkPos)) return kFALSE;
2100 TVector3 trkPosVec(trkPos[0],trkPos[1],trkPos[2]);
2101 eta = trkPosVec.Eta();
2102 phi = trkPosVec.Phi();
2103 pt = trkParam->Pt();
2105 phi += 2*TMath::Pi();
2110 //-----------------------------------------------------------------------------------
2111 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToPosition(AliExternalTrackParam *trkParam,
2112 const Float_t *clsPos,
2119 //Return the residual by extrapolating a track param to a global position
2123 if (!trkParam) return kFALSE;
2124 Double_t trkPos[3] = {0.,0.,0.};
2125 TVector3 vec(clsPos[0],clsPos[1],clsPos[2]);
2126 Double_t alpha = ((int)(vec.Phi()*TMath::RadToDeg()/20)+0.5)*20*TMath::DegToRad();
2127 vec.RotateZ(-alpha); //Rotate the cluster to the local extrapolation coordinate system
2128 if (!AliTrackerBase::PropagateTrackToBxByBz(trkParam, vec.X(), mass, step,kTRUE, 0.8, -1)) return kFALSE;
2129 if (!trkParam->GetXYZ(trkPos)) return kFALSE; //Get the extrapolated global position
2131 TVector3 clsPosVec(clsPos[0],clsPos[1],clsPos[2]);
2132 TVector3 trkPosVec(trkPos[0],trkPos[1],trkPos[2]);
2134 // track cluster matching
2135 tmpPhi = clsPosVec.DeltaPhi(trkPosVec); // tmpPhi is between -pi and pi
2136 tmpEta = clsPosVec.Eta()-trkPosVec.Eta();
2141 //----------------------------------------------------------------------------------
2142 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToCluster(AliExternalTrackParam *trkParam,
2143 const AliVCluster *cluster,
2150 //Return the residual by extrapolating a track param to a cluster
2154 if (!cluster || !trkParam)
2157 Float_t clsPos[3] = {0.,0.,0.};
2158 cluster->GetPosition(clsPos);
2160 return ExtrapolateTrackToPosition(trkParam, clsPos, mass, step, tmpEta, tmpPhi);
2163 //---------------------------------------------------------------------------------
2164 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToCluster(AliExternalTrackParam *trkParam,
2165 const AliVCluster *cluster,
2170 //Return the residual by extrapolating a track param to a clusterfStepCluster
2173 return ExtrapolateTrackToCluster(trkParam, cluster, fMass, fStepCluster, tmpEta, tmpPhi);
2176 //_______________________________________________________________________
2177 void AliEMCALRecoUtils::GetMatchedResiduals(Int_t clsIndex,
2178 Float_t &dEta, Float_t &dPhi)
2180 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2181 //Get the residuals dEta and dPhi for this cluster to the closest track
2182 //Works with ESDs and AODs
2184 if (FindMatchedPosForCluster(clsIndex) >= 999) {
2185 AliDebug(2,"No matched tracks found!\n");
2190 dEta = fResidualEta->At(FindMatchedPosForCluster(clsIndex));
2191 dPhi = fResidualPhi->At(FindMatchedPosForCluster(clsIndex));
2194 //______________________________________________________________________________________________
2195 void AliEMCALRecoUtils::GetMatchedClusterResiduals(Int_t trkIndex, Float_t &dEta, Float_t &dPhi)
2197 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2198 //Get the residuals dEta and dPhi for this track to the closest cluster
2199 //Works with ESDs and AODs
2201 if (FindMatchedPosForTrack(trkIndex) >= 999) {
2202 AliDebug(2,"No matched cluster found!\n");
2207 dEta = fResidualEta->At(FindMatchedPosForTrack(trkIndex));
2208 dPhi = fResidualPhi->At(FindMatchedPosForTrack(trkIndex));
2211 //__________________________________________________________
2212 Int_t AliEMCALRecoUtils::GetMatchedTrackIndex(Int_t clsIndex)
2214 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2215 //Get the index of matched track to this cluster
2216 //Works with ESDs and AODs
2218 if (IsClusterMatched(clsIndex))
2219 return fMatchedTrackIndex->At(FindMatchedPosForCluster(clsIndex));
2224 //__________________________________________________________
2225 Int_t AliEMCALRecoUtils::GetMatchedClusterIndex(Int_t trkIndex)
2227 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2228 //Get the index of matched cluster to this track
2229 //Works with ESDs and AODs
2231 if (IsTrackMatched(trkIndex))
2232 return fMatchedClusterIndex->At(FindMatchedPosForTrack(trkIndex));
2237 //______________________________________________________________
2238 Bool_t AliEMCALRecoUtils::IsClusterMatched(Int_t clsIndex) const
2240 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2241 //Returns if the cluster has a match
2242 if (FindMatchedPosForCluster(clsIndex) < 999)
2248 //____________________________________________________________
2249 Bool_t AliEMCALRecoUtils::IsTrackMatched(Int_t trkIndex) const
2251 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2252 //Returns if the track has a match
2253 if (FindMatchedPosForTrack(trkIndex) < 999)
2259 //______________________________________________________________________
2260 UInt_t AliEMCALRecoUtils::FindMatchedPosForCluster(Int_t clsIndex) const
2262 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2263 //Returns the position of the match in the fMatchedClusterIndex array
2264 Float_t tmpR = fCutR;
2267 for (Int_t i=0; i<fMatchedClusterIndex->GetSize(); i++)
2269 if (fMatchedClusterIndex->At(i)==clsIndex) {
2270 Float_t r = TMath::Sqrt(fResidualEta->At(i)*fResidualEta->At(i) + fResidualPhi->At(i)*fResidualPhi->At(i));
2274 AliDebug(3,Form("Matched cluster index: index: %d, dEta: %2.4f, dPhi: %2.4f.\n",
2275 fMatchedClusterIndex->At(i),fResidualEta->At(i),fResidualPhi->At(i)));
2282 //____________________________________________________________________
2283 UInt_t AliEMCALRecoUtils::FindMatchedPosForTrack(Int_t trkIndex) const
2285 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2286 //Returns the position of the match in the fMatchedTrackIndex array
2287 Float_t tmpR = fCutR;
2290 for (Int_t i=0; i<fMatchedTrackIndex->GetSize(); i++)
2292 if (fMatchedTrackIndex->At(i)==trkIndex) {
2293 Float_t r = TMath::Sqrt(fResidualEta->At(i)*fResidualEta->At(i) + fResidualPhi->At(i)*fResidualPhi->At(i));
2297 AliDebug(3,Form("Matched track index: index: %d, dEta: %2.4f, dPhi: %2.4f.\n",
2298 fMatchedTrackIndex->At(i),fResidualEta->At(i),fResidualPhi->At(i)));
2305 //__________________________________________________________________________
2306 Bool_t AliEMCALRecoUtils::IsGoodCluster(AliVCluster *cluster,
2307 const AliEMCALGeometry *geom,
2308 AliVCaloCells* cells, Int_t bc)
2310 // check if the cluster survives some quality cut
2313 Bool_t isGood=kTRUE;
2315 if (!cluster || !cluster->IsEMCAL()) return kFALSE;
2316 if (ClusterContainsBadChannel(geom,cluster->GetCellsAbsId(),cluster->GetNCells())) return kFALSE;
2317 if (!CheckCellFiducialRegion(geom,cluster,cells)) return kFALSE;
2318 if (IsExoticCluster(cluster, cells,bc)) return kFALSE;
2323 //__________________________________________________________
2324 Bool_t AliEMCALRecoUtils::IsAccepted(AliESDtrack *esdTrack)
2326 // Given a esd track, return whether the track survive all the cuts
2328 // The different quality parameter are first
2329 // retrieved from the track. then it is found out what cuts the
2330 // track did not survive and finally the cuts are imposed.
2332 UInt_t status = esdTrack->GetStatus();
2334 Int_t nClustersITS = esdTrack->GetITSclusters(0);
2335 Int_t nClustersTPC = esdTrack->GetTPCclusters(0);
2337 Float_t chi2PerClusterITS = -1;
2338 Float_t chi2PerClusterTPC = -1;
2339 if (nClustersITS!=0)
2340 chi2PerClusterITS = esdTrack->GetITSchi2()/Float_t(nClustersITS);
2341 if (nClustersTPC!=0)
2342 chi2PerClusterTPC = esdTrack->GetTPCchi2()/Float_t(nClustersTPC);
2346 if (fTrackCutsType==kGlobalCut) {
2347 Float_t maxDCAToVertexXYPtDep = 0.0182 + 0.0350/TMath::Power(esdTrack->Pt(),1.01); //This expression comes from AliESDtrackCuts::GetStandardITSTPCTrackCuts2010()
2348 //AliDebug(3,Form("Track pT = %f, DCAtoVertexXY = %f",esdTrack->Pt(),MaxDCAToVertexXYPtDep));
2349 SetMaxDCAToVertexXY(maxDCAToVertexXYPtDep); //Set pT dependent DCA cut to vertex in x-y plane
2354 esdTrack->GetImpactParameters(b,bCov);
2355 if (bCov[0]<=0 || bCov[2]<=0) {
2356 AliDebug(1, "Estimated b resolution lower or equal zero!");
2357 bCov[0]=0; bCov[2]=0;
2360 Float_t dcaToVertexXY = b[0];
2361 Float_t dcaToVertexZ = b[1];
2362 Float_t dcaToVertex = -1;
2364 if (fCutDCAToVertex2D)
2365 dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY/fCutMaxDCAToVertexXY/fCutMaxDCAToVertexXY + dcaToVertexZ*dcaToVertexZ/fCutMaxDCAToVertexZ/fCutMaxDCAToVertexZ);
2367 dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY + dcaToVertexZ*dcaToVertexZ);
2371 Bool_t cuts[kNCuts];
2372 for (Int_t i=0; i<kNCuts; i++) cuts[i]=kFALSE;
2374 // track quality cuts
2375 if (fCutRequireTPCRefit && (status&AliESDtrack::kTPCrefit)==0)
2377 if (fCutRequireITSRefit && (status&AliESDtrack::kITSrefit)==0)
2379 if (nClustersTPC<fCutMinNClusterTPC)
2381 if (nClustersITS<fCutMinNClusterITS)
2383 if (chi2PerClusterTPC>fCutMaxChi2PerClusterTPC)
2385 if (chi2PerClusterITS>fCutMaxChi2PerClusterITS)
2387 if (!fCutAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0)
2389 if (fCutDCAToVertex2D && dcaToVertex > 1)
2391 if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexXY) > fCutMaxDCAToVertexXY)
2393 if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexZ) > fCutMaxDCAToVertexZ)
2396 if (fTrackCutsType==kGlobalCut) {
2397 //Require at least one SPD point + anything else in ITS
2398 if ( (esdTrack->HasPointOnITSLayer(0) || esdTrack->HasPointOnITSLayer(1)) == kFALSE)
2403 if (fCutRequireITSStandAlone || fCutRequireITSpureSA) {
2404 if ((status & AliESDtrack::kITSin) == 0 || (status & AliESDtrack::kTPCin)) {
2408 // ITS standalone tracks
2409 if (fCutRequireITSStandAlone && !fCutRequireITSpureSA) {
2410 if (status & AliESDtrack::kITSpureSA) cuts[11] = kTRUE;
2411 } else if (fCutRequireITSpureSA) {
2412 if (!(status & AliESDtrack::kITSpureSA)) cuts[11] = kTRUE;
2418 for (Int_t i=0; i<kNCuts; i++)
2419 if (cuts[i]) { cut = kTRUE ; }
2428 //_____________________________________
2429 void AliEMCALRecoUtils::InitTrackCuts()
2431 //Intilize the track cut criteria
2432 //By default these cuts are set according to AliESDtrackCuts::GetStandardTPCOnlyTrackCuts()
2433 //Also you can customize the cuts using the setters
2435 switch (fTrackCutsType)
2439 AliInfo(Form("Track cuts for matching: GetStandardTPCOnlyTrackCuts()"));
2441 SetMinNClustersTPC(70);
2442 SetMaxChi2PerClusterTPC(4);
2443 SetAcceptKinkDaughters(kFALSE);
2444 SetRequireTPCRefit(kFALSE);
2447 SetRequireITSRefit(kFALSE);
2448 SetMaxDCAToVertexZ(3.2);
2449 SetMaxDCAToVertexXY(2.4);
2450 SetDCAToVertex2D(kTRUE);
2457 AliInfo(Form("Track cuts for matching: GetStandardITSTPCTrackCuts2010(kTURE)"));
2459 SetMinNClustersTPC(70);
2460 SetMaxChi2PerClusterTPC(4);
2461 SetAcceptKinkDaughters(kFALSE);
2462 SetRequireTPCRefit(kTRUE);
2465 SetRequireITSRefit(kTRUE);
2466 SetMaxDCAToVertexZ(2);
2467 SetMaxDCAToVertexXY();
2468 SetDCAToVertex2D(kFALSE);
2475 AliInfo(Form("Track cuts for matching: Loose cut w/o DCA cut"));
2476 SetMinNClustersTPC(50);
2477 SetAcceptKinkDaughters(kTRUE);
2482 case kITSStandAlone:
2484 AliInfo(Form("Track cuts for matching: ITS Stand Alone tracks cut w/o DCA cut"));
2485 SetRequireITSRefit(kTRUE);
2486 SetRequireITSStandAlone(kTRUE);
2487 SetITSTrackSA(kTRUE);
2495 //________________________________________________________________________
2496 void AliEMCALRecoUtils::SetClusterMatchedToTrack(const AliVEvent *event)
2498 // Checks if tracks are matched to EMC clusters and set the matched EMCAL cluster index to ESD track.
2500 Int_t nTracks = event->GetNumberOfTracks();
2501 for (Int_t iTrack = 0; iTrack < nTracks; ++iTrack)
2503 AliVTrack* track = dynamic_cast<AliVTrack*>(event->GetTrack(iTrack));
2506 AliWarning(Form("Could not receive track %d", iTrack));
2510 Int_t matchClusIndex = GetMatchedClusterIndex(iTrack);
2511 track->SetEMCALcluster(matchClusIndex); //sets -1 if track not matched within residual
2512 /*the following can be done better if AliVTrack::SetStatus will be there. Patch pending with Andreas/Peter*/
2513 AliESDtrack* esdtrack = dynamic_cast<AliESDtrack*>(track);
2515 if (matchClusIndex != -1)
2516 esdtrack->SetStatus(AliESDtrack::kEMCALmatch);
2518 esdtrack->ResetStatus(AliESDtrack::kEMCALmatch);
2520 AliAODTrack* aodtrack = dynamic_cast<AliAODTrack*>(track);
2521 if (matchClusIndex != -1)
2522 aodtrack->SetStatus(AliESDtrack::kEMCALmatch);
2524 aodtrack->ResetStatus(AliESDtrack::kEMCALmatch);
2527 AliDebug(2,"Track matched to closest cluster");
2530 //_________________________________________________________________________
2531 void AliEMCALRecoUtils::SetTracksMatchedToCluster(const AliVEvent *event)
2533 // Checks if EMC clusters are matched to ESD track.
2534 // Adds track indexes of all the tracks matched to a cluster withing residuals in ESDCalocluster.
2536 for (Int_t iClus=0; iClus < event->GetNumberOfCaloClusters(); ++iClus)
2538 AliVCluster *cluster = event->GetCaloCluster(iClus);
2539 if (!cluster->IsEMCAL())
2542 Int_t nTracks = event->GetNumberOfTracks();
2543 TArrayI arrayTrackMatched(nTracks);
2545 // Get the closest track matched to the cluster
2547 Int_t matchTrackIndex = GetMatchedTrackIndex(iClus);
2548 if (matchTrackIndex != -1)
2550 arrayTrackMatched[nMatched] = matchTrackIndex;
2554 // Get all other tracks matched to the cluster
2555 for (Int_t iTrk=0; iTrk<nTracks; ++iTrk)
2557 AliVTrack* track = dynamic_cast<AliVTrack*>(event->GetTrack(iTrk));
2559 if( !track ) continue;
2561 if ( iTrk == matchTrackIndex ) continue;
2563 if ( track->GetEMCALcluster() == iClus )
2565 arrayTrackMatched[nMatched] = iTrk;
2570 //printf("Tender::SetTracksMatchedToCluster - cluster E %f, N matches %d, first match %d\n",cluster->E(),nMatched,arrayTrackMatched[0]);
2572 arrayTrackMatched.Set(nMatched);
2573 AliESDCaloCluster *esdcluster = dynamic_cast<AliESDCaloCluster*>(cluster);
2575 esdcluster->AddTracksMatched(arrayTrackMatched);
2576 else if (nMatched>0) {
2577 AliAODCaloCluster *aodcluster = dynamic_cast<AliAODCaloCluster*>(cluster);
2579 aodcluster->AddTrackMatched(event->GetTrack(arrayTrackMatched.At(0)));
2582 Float_t eta= -999, phi = -999;
2583 if (matchTrackIndex != -1)
2584 GetMatchedResiduals(iClus, eta, phi);
2585 cluster->SetTrackDistance(phi, eta);
2588 AliDebug(2,"Cluster matched to tracks");
2591 //___________________________________________________
2592 void AliEMCALRecoUtils::Print(const Option_t *) const
2596 printf("AliEMCALRecoUtils Settings: \n");
2597 printf("Misalignment shifts\n");
2598 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,
2599 fMisalTransShift[i*3],fMisalTransShift[i*3+1],fMisalTransShift[i*3+2],
2600 fMisalRotShift[i*3], fMisalRotShift[i*3+1], fMisalRotShift[i*3+2] );
2601 printf("Non linearity function %d, parameters:\n", fNonLinearityFunction);
2602 for (Int_t i=0; i<6; i++) printf("param[%d]=%f\n",i, fNonLinearityParams[i]);
2604 printf("Position Recalculation option %d, Particle Type %d, fW0 %2.2f, Recalibrate Data %d \n",fPosAlgo,fParticleType,fW0, fRecalibration);
2606 printf("Matching criteria: ");
2607 if (fCutEtaPhiSum) {
2608 printf("sqrt(dEta^2+dPhi^2)<%4.3f\n",fCutR);
2609 } else if (fCutEtaPhiSeparate) {
2610 printf("dEta<%4.3f, dPhi<%4.3f\n",fCutEta,fCutPhi);
2613 printf("please specify your cut criteria\n");
2614 printf("To cut on sqrt(dEta^2+dPhi^2), use: SwitchOnCutEtaPhiSum()\n");
2615 printf("To cut on dEta and dPhi separately, use: SwitchOnCutEtaPhiSeparate()\n");
2618 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);
2619 printf("Cluster selection window: dR < %2.0f\n",fClusterWindow);
2621 printf("Track cuts: \n");
2622 printf("Minimum track pT: %1.2f\n",fCutMinTrackPt);
2623 printf("AOD track selection: tpc only %d, or hybrid %d, or mask: %d\n",fAODTPCOnlyTracks,fAODHybridTracks, fAODFilterMask);
2624 printf("TPCRefit = %d, ITSRefit = %d\n",fCutRequireTPCRefit,fCutRequireITSRefit);
2625 printf("AcceptKinks = %d\n",fCutAcceptKinkDaughters);
2626 printf("MinNCulsterTPC = %d, MinNClusterITS = %d\n",fCutMinNClusterTPC,fCutMinNClusterITS);
2627 printf("MaxChi2TPC = %2.2f, MaxChi2ITS = %2.2f\n",fCutMaxChi2PerClusterTPC,fCutMaxChi2PerClusterITS);
2628 printf("DCSToVertex2D = %d, MaxDCAToVertexXY = %2.2f, MaxDCAToVertexZ = %2.2f\n",fCutDCAToVertex2D,fCutMaxDCAToVertexXY,fCutMaxDCAToVertexZ);