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] ; }
151 //______________________________________________________________________
152 AliEMCALRecoUtils & AliEMCALRecoUtils::operator = (const AliEMCALRecoUtils & reco)
154 //Assignment operator
156 if (this == &reco)return *this;
157 ((TNamed *)this)->operator=(reco);
159 for (Int_t i = 0; i < 15 ; i++) { fMisalTransShift[i] = reco.fMisalTransShift[i] ;
160 fMisalRotShift[i] = reco.fMisalRotShift[i] ; }
161 for (Int_t i = 0; i < 7 ; i++) { fNonLinearityParams[i] = reco.fNonLinearityParams[i] ; }
162 for (Int_t i = 0; i < 3 ; i++) { fSmearClusterParam[i] = reco.fSmearClusterParam[i] ; }
164 fParticleType = reco.fParticleType;
165 fPosAlgo = reco.fPosAlgo;
168 fNonLinearityFunction = reco.fNonLinearityFunction;
169 fNonLinearThreshold = reco.fNonLinearThreshold;
170 fSmearClusterEnergy = reco.fSmearClusterEnergy;
172 fCellsRecalibrated = reco.fCellsRecalibrated;
173 fRecalibration = reco.fRecalibration;
174 fEMCALRecalibrationFactors = reco.fEMCALRecalibrationFactors;
176 fTimeRecalibration = reco.fTimeRecalibration;
177 fEMCALTimeRecalibrationFactors = reco.fEMCALTimeRecalibrationFactors;
179 fUseRunCorrectionFactors = reco.fUseRunCorrectionFactors;
181 fRemoveBadChannels = reco.fRemoveBadChannels;
182 fRecalDistToBadChannels = reco.fRecalDistToBadChannels;
183 fEMCALBadChannelMap = reco.fEMCALBadChannelMap;
185 fNCellsFromEMCALBorder = reco.fNCellsFromEMCALBorder;
186 fNoEMCALBorderAtEta0 = reco.fNoEMCALBorderAtEta0;
188 fRejectExoticCluster = reco.fRejectExoticCluster;
189 fRejectExoticCells = reco.fRejectExoticCells;
190 fExoticCellFraction = reco.fExoticCellFraction;
191 fExoticCellDiffTime = reco.fExoticCellDiffTime;
192 fExoticCellMinAmplitude = reco.fExoticCellMinAmplitude;
194 fPIDUtils = reco.fPIDUtils;
196 fAODFilterMask = reco.fAODFilterMask;
197 fAODHybridTracks = reco.fAODHybridTracks;
198 fAODTPCOnlyTracks = reco.fAODTPCOnlyTracks;
200 fCutEtaPhiSum = reco.fCutEtaPhiSum;
201 fCutEtaPhiSeparate = reco.fCutEtaPhiSeparate;
203 fCutEta = reco.fCutEta;
204 fCutPhi = reco.fCutPhi;
205 fClusterWindow = reco.fClusterWindow;
207 fStepSurface = reco.fStepSurface;
208 fStepCluster = reco.fStepCluster;
209 fITSTrackSA = reco.fITSTrackSA;
210 fEMCalSurfaceDistance = reco.fEMCalSurfaceDistance;
212 fTrackCutsType = reco.fTrackCutsType;
213 fCutMinTrackPt = reco.fCutMinTrackPt;
214 fCutMinNClusterTPC = reco.fCutMinNClusterTPC;
215 fCutMinNClusterITS = reco.fCutMinNClusterITS;
216 fCutMaxChi2PerClusterTPC = reco.fCutMaxChi2PerClusterTPC;
217 fCutMaxChi2PerClusterITS = reco.fCutMaxChi2PerClusterITS;
218 fCutRequireTPCRefit = reco.fCutRequireTPCRefit;
219 fCutRequireITSRefit = reco.fCutRequireITSRefit;
220 fCutAcceptKinkDaughters = reco.fCutAcceptKinkDaughters;
221 fCutMaxDCAToVertexXY = reco.fCutMaxDCAToVertexXY;
222 fCutMaxDCAToVertexZ = reco.fCutMaxDCAToVertexZ;
223 fCutDCAToVertex2D = reco.fCutDCAToVertex2D;
224 fCutRequireITSStandAlone = reco.fCutRequireITSStandAlone;
225 fCutRequireITSpureSA = reco.fCutRequireITSpureSA;
227 if (reco.fResidualEta) {
228 // assign or copy construct
230 *fResidualEta = *reco.fResidualEta;
232 fResidualEta = new TArrayF(*reco.fResidualEta);
235 if (fResidualEta) delete fResidualEta;
239 if (reco.fResidualPhi) {
240 // assign or copy construct
242 *fResidualPhi = *reco.fResidualPhi;
244 fResidualPhi = new TArrayF(*reco.fResidualPhi);
247 if (fResidualPhi) delete fResidualPhi;
251 if (reco.fMatchedTrackIndex) {
252 // assign or copy construct
253 if (fMatchedTrackIndex) {
254 *fMatchedTrackIndex = *reco.fMatchedTrackIndex;
256 fMatchedTrackIndex = new TArrayI(*reco.fMatchedTrackIndex);
259 if (fMatchedTrackIndex) delete fMatchedTrackIndex;
260 fMatchedTrackIndex = 0;
263 if (reco.fMatchedClusterIndex) {
264 // assign or copy construct
265 if (fMatchedClusterIndex) {
266 *fMatchedClusterIndex = *reco.fMatchedClusterIndex;
268 fMatchedClusterIndex = new TArrayI(*reco.fMatchedClusterIndex);
271 if (fMatchedClusterIndex) delete fMatchedClusterIndex;
272 fMatchedClusterIndex = 0;
278 //_____________________________________
279 AliEMCALRecoUtils::~AliEMCALRecoUtils()
283 if (fEMCALRecalibrationFactors) {
284 fEMCALRecalibrationFactors->Clear();
285 delete fEMCALRecalibrationFactors;
288 if (fEMCALTimeRecalibrationFactors) {
289 fEMCALTimeRecalibrationFactors->Clear();
290 delete fEMCALTimeRecalibrationFactors;
293 if (fEMCALBadChannelMap) {
294 fEMCALBadChannelMap->Clear();
295 delete fEMCALBadChannelMap;
298 delete fMatchedTrackIndex ;
299 delete fMatchedClusterIndex ;
300 delete fResidualEta ;
301 delete fResidualPhi ;
307 //_______________________________________________________________________________
308 Bool_t AliEMCALRecoUtils::AcceptCalibrateCell(Int_t absID, Int_t bc,
309 Float_t & amp, Double_t & time,
310 AliVCaloCells* cells)
312 // Reject cell if criteria not passed and calibrate it
314 AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
316 if (absID < 0 || absID >= 24*48*geom->GetNumberOfSuperModules())
319 Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1;
321 if (!geom->GetCellIndex(absID,imod,iTower,iIphi,iIeta)) {
322 // cell absID does not exist
327 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);
329 // Do not include bad channels found in analysis,
330 if (IsBadChannelsRemovalSwitchedOn() && GetEMCALChannelStatus(imod, ieta, iphi)) {
335 amp = cells->GetCellAmplitude(absID);
336 if (!fCellsRecalibrated && IsRecalibrationOn())
337 amp *= GetEMCALChannelRecalibrationFactor(imod,ieta,iphi);
340 time = cells->GetCellTime(absID);
342 RecalibrateCellTime(absID,bc,time);
347 //_____________________________________________________________________________
348 Bool_t AliEMCALRecoUtils::CheckCellFiducialRegion(const AliEMCALGeometry* geom,
349 const AliVCluster* cluster,
350 AliVCaloCells* cells)
352 // Given the list of AbsId of the cluster, get the maximum cell and
353 // check if there are fNCellsFromBorder from the calorimeter border
357 AliInfo("Cluster pointer null!");
361 //If the distance to the border is 0 or negative just exit accept all clusters
362 if (cells->GetType()==AliVCaloCells::kEMCALCell && fNCellsFromEMCALBorder <= 0 )
365 Int_t absIdMax = -1, iSM =-1, ieta = -1, iphi = -1;
366 Bool_t shared = kFALSE;
367 GetMaxEnergyCell(geom, cells, cluster, absIdMax, iSM, ieta, iphi, shared);
369 AliDebug(2,Form("Cluster Max AbsId %d, Cell Energy %2.2f, Cluster Energy %2.2f, Ncells from border %d, EMCAL eta=0 %d\n",
370 absIdMax, cells->GetCellAmplitude(absIdMax), cluster->E(), fNCellsFromEMCALBorder, fNoEMCALBorderAtEta0));
372 if (absIdMax==-1) return kFALSE;
374 //Check if the cell is close to the borders:
375 Bool_t okrow = kFALSE;
376 Bool_t okcol = kFALSE;
378 if (iSM < 0 || iphi < 0 || ieta < 0 ) {
379 AliFatal(Form("Negative value for super module: %d, or cell ieta: %d, or cell iphi: %d, check EMCAL geometry name\n",
381 return kFALSE; // trick coverity
386 if( geom->GetSMType(iSM) == AliEMCALGeometry::kEMCAL_Half ) iPhiLast /= 2;
387 else if ( geom->GetSMType(iSM) == AliEMCALGeometry::kEMCAL_3rd ) iPhiLast /= 3;// 1/3 sm case
389 if(iphi >= fNCellsFromEMCALBorder && iphi < iPhiLast - fNCellsFromEMCALBorder) okrow = kTRUE;
393 if(!fNoEMCALBorderAtEta0 || geom->IsDCALSM(iSM)) {// conside inner border
394 if( geom->GetSMType(iSM) == AliEMCALGeometry::kDCAL_Standard ) iEtaLast = iEtaLast*2/3;
395 if(ieta > fNCellsFromEMCALBorder && ieta < iEtaLast-fNCellsFromEMCALBorder) okcol = kTRUE;
398 if (ieta >= fNCellsFromEMCALBorder) okcol = kTRUE;
400 if(ieta < iEtaLast-fNCellsFromEMCALBorder) okcol = kTRUE;
404 AliDebug(2,Form("EMCAL Cluster in %d cells fiducial volume: ieta %d, iphi %d, SM %d: column? %d, row? %d\nq",
405 fNCellsFromEMCALBorder, ieta, iphi, iSM, okcol, okrow));
407 if (okcol && okrow) {
408 //printf("Accept\n");
411 //printf("Reject\n");
412 AliDebug(2,Form("Reject cluster in border, max cell : ieta %d, iphi %d, SM %d\n",ieta, iphi, iSM));
417 //_______________________________________________________________________________
418 Bool_t AliEMCALRecoUtils::ClusterContainsBadChannel(const AliEMCALGeometry* geom,
419 const UShort_t* cellList,
422 // Check that in the cluster cells, there is no bad channel of those stored
423 // in fEMCALBadChannelMap or fPHOSBadChannelMap
425 if (!fRemoveBadChannels) return kFALSE;
426 if (!fEMCALBadChannelMap) return kFALSE;
431 for (Int_t iCell = 0; iCell<nCells; iCell++) {
432 //Get the column and row
433 Int_t iTower = -1, iIphi = -1, iIeta = -1;
434 geom->GetCellIndex(cellList[iCell],imod,iTower,iIphi,iIeta);
435 if (fEMCALBadChannelMap->GetEntries() <= imod) continue;
436 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,irow,icol);
437 if (GetEMCALChannelStatus(imod, icol, irow)) {
438 AliDebug(2,Form("Cluster with bad channel: SM %d, col %d, row %d\n",imod, icol, irow));
441 }// cell cluster loop
447 //___________________________________________________________________________
448 Float_t AliEMCALRecoUtils::GetECross(Int_t absID, Double_t tcell,
449 AliVCaloCells* cells, Int_t bc)
451 //Calculate the energy in the cross around the energy given cell
453 AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
455 Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1;
456 geom->GetCellIndex(absID,imod,iTower,iIphi,iIeta);
457 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);
459 //Get close cells index, energy and time, not in corners
464 if ( iphi < AliEMCALGeoParams::fgkEMCALRows-1) absID1 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi+1, ieta);
465 if ( iphi > 0 ) absID2 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi-1, ieta);
467 // In case of cell in eta = 0 border, depending on SM shift the cross cell index
472 if ( ieta == AliEMCALGeoParams::fgkEMCALCols-1 && !(imod%2) ) {
473 absID3 = geom-> GetAbsCellIdFromCellIndexes(imod+1, iphi, 0);
474 absID4 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta-1);
475 } else if ( ieta == 0 && imod%2 ) {
476 absID3 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta+1);
477 absID4 = geom-> GetAbsCellIdFromCellIndexes(imod-1, iphi, AliEMCALGeoParams::fgkEMCALCols-1);
479 if ( ieta < AliEMCALGeoParams::fgkEMCALCols-1 )
480 absID3 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta+1);
482 absID4 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta-1);
485 //printf("IMOD %d, AbsId %d, a %d, b %d, c %d e %d \n",imod,absID,absID1,absID2,absID3,absID4);
487 Float_t ecell1 = 0, ecell2 = 0, ecell3 = 0, ecell4 = 0;
488 Double_t tcell1 = 0, tcell2 = 0, tcell3 = 0, tcell4 = 0;
490 AcceptCalibrateCell(absID1,bc, ecell1,tcell1,cells);
491 AcceptCalibrateCell(absID2,bc, ecell2,tcell2,cells);
492 AcceptCalibrateCell(absID3,bc, ecell3,tcell3,cells);
493 AcceptCalibrateCell(absID4,bc, ecell4,tcell4,cells);
495 if (TMath::Abs(tcell-tcell1)*1.e9 > fExoticCellDiffTime) ecell1 = 0 ;
496 if (TMath::Abs(tcell-tcell2)*1.e9 > fExoticCellDiffTime) ecell2 = 0 ;
497 if (TMath::Abs(tcell-tcell3)*1.e9 > fExoticCellDiffTime) ecell3 = 0 ;
498 if (TMath::Abs(tcell-tcell4)*1.e9 > fExoticCellDiffTime) ecell4 = 0 ;
500 return ecell1+ecell2+ecell3+ecell4;
503 //_____________________________________________________________________________________________
504 Bool_t AliEMCALRecoUtils::IsExoticCell(Int_t absID, AliVCaloCells* cells, Int_t bc)
506 // Look to cell neighbourhood and reject if it seems exotic
507 // Do before recalibrating the cells
509 if (!fRejectExoticCells) return kFALSE;
513 Bool_t accept = AcceptCalibrateCell(absID, bc, ecell ,tcell ,cells);
515 if (!accept) return kTRUE; // reject this cell
517 if (ecell < fExoticCellMinAmplitude) return kFALSE; // do not reject low energy cells
519 Float_t eCross = GetECross(absID,tcell,cells,bc);
521 if (1-eCross/ecell > fExoticCellFraction) {
522 AliDebug(2,Form("AliEMCALRecoUtils::IsExoticCell() - EXOTIC CELL id %d, eCell %f, eCross %f, 1-eCross/eCell %f\n",
523 absID,ecell,eCross,1-eCross/ecell));
530 //___________________________________________________________________
531 Bool_t AliEMCALRecoUtils::IsExoticCluster(const AliVCluster *cluster,
532 AliVCaloCells *cells,
535 // Check if the cluster highest energy tower is exotic
538 AliInfo("Cluster pointer null!");
542 if (!fRejectExoticCluster) return kFALSE;
544 // Get highest energy tower
545 AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
546 Int_t iSupMod = -1, absId = -1, ieta = -1, iphi = -1;
547 Bool_t shared = kFALSE;
548 GetMaxEnergyCell(geom, cells, cluster, absId, iSupMod, ieta, iphi, shared);
550 return IsExoticCell(absId,cells,bc);
553 //_______________________________________________________________________
554 Float_t AliEMCALRecoUtils::SmearClusterEnergy(const AliVCluster* cluster)
556 //In case of MC analysis, smear energy to match resolution/calibration in real data
559 AliInfo("Cluster pointer null!");
563 Float_t energy = cluster->E() ;
564 Float_t rdmEnergy = energy ;
565 if (fSmearClusterEnergy) {
566 rdmEnergy = fRandom.Gaus(energy,fSmearClusterParam[0] * TMath::Sqrt(energy) +
567 fSmearClusterParam[1] * energy +
568 fSmearClusterParam[2] );
569 AliDebug(2, Form("Energy: original %f, smeared %f\n", energy, rdmEnergy));
575 //____________________________________________________________________________
576 Float_t AliEMCALRecoUtils::CorrectClusterEnergyLinearity(AliVCluster* cluster)
578 // Correct cluster energy from non linearity functions
581 AliInfo("Cluster pointer null!");
585 Float_t energy = cluster->E();
588 // Clusters with less than 50 MeV or negative are not possible
589 AliInfo(Form("Too Low Cluster energy!, E = %f < 0.05 GeV",energy));
593 switch (fNonLinearityFunction)
597 //Non-Linearity correction (from MC with function ([0]*exp(-[1]/E))+(([2]/([3]*2.*TMath::Pi())*exp(-(E-[4])^2/(2.*[3]^2)))))
598 //fNonLinearityParams[0] = 1.014;
599 //fNonLinearityParams[1] =-0.03329;
600 //fNonLinearityParams[2] =-0.3853;
601 //fNonLinearityParams[3] = 0.5423;
602 //fNonLinearityParams[4] =-0.4335;
603 energy *= (fNonLinearityParams[0]*exp(-fNonLinearityParams[1]/energy))+
604 ((fNonLinearityParams[2]/(fNonLinearityParams[3]*2.*TMath::Pi())*
605 exp(-(energy-fNonLinearityParams[4])*(energy-fNonLinearityParams[4])/(2.*fNonLinearityParams[3]*fNonLinearityParams[3]))));
611 //Non-Linearity correction (from MC with function [0]/((x+[1])^[2]))+1;
612 //fNonLinearityParams[0] = 3.11111e-02;
613 //fNonLinearityParams[1] =-5.71666e-02;
614 //fNonLinearityParams[2] = 5.67995e-01;
616 energy *= fNonLinearityParams[0]/TMath::Power(energy+fNonLinearityParams[1],fNonLinearityParams[2])+1;
622 //Same as beam test corrected, change parameters
623 //fNonLinearityParams[0] = 9.81039e-01
624 //fNonLinearityParams[1] = 1.13508e-01;
625 //fNonLinearityParams[2] = 1.00173e+00;
626 //fNonLinearityParams[3] = 9.67998e-02;
627 //fNonLinearityParams[4] = 2.19381e+02;
628 //fNonLinearityParams[5] = 6.31604e+01;
629 //fNonLinearityParams[6] = 1;
630 energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
638 //Non-Linearity correction (from Olga Data with function p0+p1*exp(-p2*E))
639 //fNonLinearityParams[0] = 1.04;
640 //fNonLinearityParams[1] = -0.1445;
641 //fNonLinearityParams[2] = 1.046;
642 energy /= (fNonLinearityParams[0]+fNonLinearityParams[1]*exp(-fNonLinearityParams[2]*energy)); //Olga function
646 case kPi0GammaConversion:
648 //Non-Linearity correction (Nicolas from Dimitri Data with function C*[1-a*exp(-b*E)])
649 //fNonLinearityParams[0] = 0.139393/0.1349766;
650 //fNonLinearityParams[1] = 0.0566186;
651 //fNonLinearityParams[2] = 0.982133;
652 energy /= fNonLinearityParams[0]*(1-fNonLinearityParams[1]*exp(-fNonLinearityParams[2]*energy));
659 //From beam test, Alexei's results, for different ZS thresholds
660 // th=30 MeV; th = 45 MeV; th = 75 MeV
661 //fNonLinearityParams[0] = 1.007; 1.003; 1.002
662 //fNonLinearityParams[1] = 0.894; 0.719; 0.797
663 //fNonLinearityParams[2] = 0.246; 0.334; 0.358
664 //Rescale the param[0] with 1.03
665 energy /= fNonLinearityParams[0]/(1+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]));
670 case kBeamTestCorrected:
672 //From beam test, corrected for material between beam and EMCAL
673 //fNonLinearityParams[0] = 0.99078
674 //fNonLinearityParams[1] = 0.161499;
675 //fNonLinearityParams[2] = 0.655166;
676 //fNonLinearityParams[3] = 0.134101;
677 //fNonLinearityParams[4] = 163.282;
678 //fNonLinearityParams[5] = 23.6904;
679 //fNonLinearityParams[6] = 0.978;
680 energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
685 case kBeamTestCorrectedv2:
687 //From beam test, corrected for material between beam and EMCAL
688 //fNonLinearityParams[0] = 0.983504;
689 //fNonLinearityParams[1] = 0.210106;
690 //fNonLinearityParams[2] = 0.897274;
691 //fNonLinearityParams[3] = 0.0829064;
692 //fNonLinearityParams[4] = 152.299;
693 //fNonLinearityParams[5] = 31.5028;
694 //fNonLinearityParams[6] = 0.968;
695 energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
702 //Based on fit to the MC/data using kNoCorrection on the data - utilizes symmetric decay method and kPi0MCv5(MC) - 28 Oct 2013
703 //fNonLinearityParams[0] = 1.0;
704 //fNonLinearityParams[1] = 6.64778e-02;
705 //fNonLinearityParams[2] = 1.570;
706 //fNonLinearityParams[3] = 9.67998e-02;
707 //fNonLinearityParams[4] = 2.19381e+02;
708 //fNonLinearityParams[5] = 6.31604e+01;
709 //fNonLinearityParams[6] = 1.01286;
710 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));
717 //Based on comparing MC truth information to the reconstructed energy of clusters.
718 //fNonLinearityParams[0] = 1.0;
719 //fNonLinearityParams[1] = 6.64778e-02;
720 //fNonLinearityParams[2] = 1.570;
721 //fNonLinearityParams[3] = 9.67998e-02;
722 //fNonLinearityParams[4] = 2.19381e+02;
723 //fNonLinearityParams[5] = 6.31604e+01;
724 //fNonLinearityParams[6] = 1.01286;
725 energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
732 //Based on fit to the MC/data using kNoCorrection on the data
733 // - utilizes symmetric decay method and kPi0MCv6(MC) - 09 Dec 2014
734 // - parameters constrained by the test beam data as well
735 // described in the note: https://aliceinfo.cern.ch/Notes/node/211 - Sec 3.1.2 (Test Beam Constrained SDM).
736 //fNonLinearityParams[0] = 1.0;
737 //fNonLinearityParams[1] = 0.237767;
738 //fNonLinearityParams[2] = 0.651203;
739 //fNonLinearityParams[3] = 0.183741;
740 //fNonLinearityParams[4] = 155.427;
741 //fNonLinearityParams[5] = 17.0335;
742 //fNonLinearityParams[6] = 0.987054;
743 energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
750 //Based on comparing MC truth information to the reconstructed energy of clusters.
751 // described in the note: https://aliceinfo.cern.ch/Notes/node/211 - Sec 3.1.2 (Test Beam Constrained SDM).
752 //fNonLinearityParams[0] = 1.0;
753 //fNonLinearityParams[1] = 0.0797873;
754 //fNonLinearityParams[2] = 1.68322;
755 //fNonLinearityParams[3] = 0.0806098;
756 //fNonLinearityParams[4] = 244.586;
757 //fNonLinearityParams[5] = 116.938;
758 //fNonLinearityParams[6] = 1.00437;
759 energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
765 AliDebug(2,"No correction on the energy\n");
773 //__________________________________________________
774 void AliEMCALRecoUtils::InitNonLinearityParam()
776 //Initialising Non Linearity Parameters
778 if (fNonLinearityFunction == kPi0MC) {
779 fNonLinearityParams[0] = 1.014;
780 fNonLinearityParams[1] = -0.03329;
781 fNonLinearityParams[2] = -0.3853;
782 fNonLinearityParams[3] = 0.5423;
783 fNonLinearityParams[4] = -0.4335;
786 if (fNonLinearityFunction == kPi0MCv2) {
787 fNonLinearityParams[0] = 3.11111e-02;
788 fNonLinearityParams[1] =-5.71666e-02;
789 fNonLinearityParams[2] = 5.67995e-01;
792 if (fNonLinearityFunction == kPi0MCv3) {
793 fNonLinearityParams[0] = 9.81039e-01;
794 fNonLinearityParams[1] = 1.13508e-01;
795 fNonLinearityParams[2] = 1.00173e+00;
796 fNonLinearityParams[3] = 9.67998e-02;
797 fNonLinearityParams[4] = 2.19381e+02;
798 fNonLinearityParams[5] = 6.31604e+01;
799 fNonLinearityParams[6] = 1;
802 if (fNonLinearityFunction == kPi0GammaGamma) {
803 fNonLinearityParams[0] = 1.04;
804 fNonLinearityParams[1] = -0.1445;
805 fNonLinearityParams[2] = 1.046;
808 if (fNonLinearityFunction == kPi0GammaConversion) {
809 fNonLinearityParams[0] = 0.139393;
810 fNonLinearityParams[1] = 0.0566186;
811 fNonLinearityParams[2] = 0.982133;
814 if (fNonLinearityFunction == kBeamTest) {
815 if (fNonLinearThreshold == 30) {
816 fNonLinearityParams[0] = 1.007;
817 fNonLinearityParams[1] = 0.894;
818 fNonLinearityParams[2] = 0.246;
820 if (fNonLinearThreshold == 45) {
821 fNonLinearityParams[0] = 1.003;
822 fNonLinearityParams[1] = 0.719;
823 fNonLinearityParams[2] = 0.334;
825 if (fNonLinearThreshold == 75) {
826 fNonLinearityParams[0] = 1.002;
827 fNonLinearityParams[1] = 0.797;
828 fNonLinearityParams[2] = 0.358;
832 if (fNonLinearityFunction == kBeamTestCorrected) {
833 fNonLinearityParams[0] = 0.99078;
834 fNonLinearityParams[1] = 0.161499;
835 fNonLinearityParams[2] = 0.655166;
836 fNonLinearityParams[3] = 0.134101;
837 fNonLinearityParams[4] = 163.282;
838 fNonLinearityParams[5] = 23.6904;
839 fNonLinearityParams[6] = 0.978;
842 if (fNonLinearityFunction == kBeamTestCorrectedv2) {
843 fNonLinearityParams[0] = 0.983504;
844 fNonLinearityParams[1] = 0.210106;
845 fNonLinearityParams[2] = 0.897274;
846 fNonLinearityParams[3] = 0.0829064;
847 fNonLinearityParams[4] = 152.299;
848 fNonLinearityParams[5] = 31.5028;
849 fNonLinearityParams[6] = 0.968;
852 if (fNonLinearityFunction == kSDMv5) {
853 fNonLinearityParams[0] = 1.0;
854 fNonLinearityParams[1] = 6.64778e-02;
855 fNonLinearityParams[2] = 1.570;
856 fNonLinearityParams[3] = 9.67998e-02;
857 fNonLinearityParams[4] = 2.19381e+02;
858 fNonLinearityParams[5] = 6.31604e+01;
859 fNonLinearityParams[6] = 1.01286;
862 if (fNonLinearityFunction == kPi0MCv5) {
863 fNonLinearityParams[0] = 1.0;
864 fNonLinearityParams[1] = 6.64778e-02;
865 fNonLinearityParams[2] = 1.570;
866 fNonLinearityParams[3] = 9.67998e-02;
867 fNonLinearityParams[4] = 2.19381e+02;
868 fNonLinearityParams[5] = 6.31604e+01;
869 fNonLinearityParams[6] = 1.01286;
872 if (fNonLinearityFunction == kSDMv6) {
873 fNonLinearityParams[0] = 1.0;
874 fNonLinearityParams[1] = 0.237767;
875 fNonLinearityParams[2] = 0.651203;
876 fNonLinearityParams[3] = 0.183741;
877 fNonLinearityParams[4] = 155.427;
878 fNonLinearityParams[5] = 17.0335;
879 fNonLinearityParams[6] = 0.987054;
882 if (fNonLinearityFunction == kPi0MCv6) {
883 fNonLinearityParams[0] = 1.0;
884 fNonLinearityParams[1] = 0.0797873;
885 fNonLinearityParams[2] = 1.68322;
886 fNonLinearityParams[3] = 0.0806098;
887 fNonLinearityParams[4] = 244.586;
888 fNonLinearityParams[5] = 116.938;
889 fNonLinearityParams[6] = 1.00437;
894 //_________________________________________________________
895 Float_t AliEMCALRecoUtils::GetDepth(Float_t energy,
899 //Calculate shower depth for a given cluster energy and particle type
905 Float_t arg = energy*1000/ ecr; //Multiply energy by 1000 to transform to MeV
913 depth = x0 * (TMath::Log(arg) + 0.5);
920 depth = x0 * (TMath::Log(arg) - 0.5);
927 gGeoManager->cd("ALIC_1/XEN1_1");
928 TGeoNode *geoXEn1 = gGeoManager->GetCurrentNode();
929 TGeoNodeMatrix *geoSM = dynamic_cast<TGeoNodeMatrix *>(geoXEn1->GetDaughter(iSM));
931 TGeoVolume *geoSMVol = geoSM->GetVolume();
932 TGeoShape *geoSMShape = geoSMVol->GetShape();
933 TGeoBBox *geoBox = dynamic_cast<TGeoBBox *>(geoSMShape);
934 if (geoBox) depth = 0.5 * geoBox->GetDX()*2 ;
935 else AliFatal("Null GEANT box");
937 else AliFatal("NULL GEANT node matrix");
944 depth = x0 * (TMath::Log(arg) - 0.5);
953 depth = x0 * (TMath::Log(arg) + 0.5);
959 //____________________________________________________________________
960 void AliEMCALRecoUtils::GetMaxEnergyCell(const AliEMCALGeometry *geom,
961 AliVCaloCells* cells,
962 const AliVCluster* clu,
969 //For a given CaloCluster gets the absId of the cell
970 //with maximum energy deposit.
973 Double_t eCell = -1.;
974 Float_t fraction = 1.;
975 Float_t recalFactor = 1.;
976 Int_t cellAbsId = -1 ;
984 AliInfo("Cluster pointer null!");
985 absId=-1; iSupMod0=-1, ieta = -1; iphi = -1; shared = -1;
989 for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++) {
990 cellAbsId = clu->GetCellAbsId(iDig);
991 fraction = clu->GetCellAmplitudeFraction(iDig);
992 //printf("a Cell %d, id, %d, amp %f, fraction %f\n",iDig,cellAbsId,cells->GetCellAmplitude(cellAbsId),fraction);
993 if (fraction < 1e-4) fraction = 1.; // in case unfolding is off
994 geom->GetCellIndex(cellAbsId,iSupMod,iTower,iIphi,iIeta);
995 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);
998 } else if (iSupMod0!=iSupMod) {
1000 //printf("AliEMCALRecoUtils::GetMaxEnergyCell() - SHARED CLUSTER\n");
1002 if (!fCellsRecalibrated && IsRecalibrationOn()) {
1003 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1005 eCell = cells->GetCellAmplitude(cellAbsId)*fraction*recalFactor;
1006 //printf("b Cell %d, id, %d, amp %f, fraction %f\n",iDig,cellAbsId,eCell,fraction);
1010 //printf("\t new max: cell %d, e %f, ecell %f\n",maxId, eMax,eCell);
1014 //Get from the absid the supermodule, tower and eta/phi numbers
1015 geom->GetCellIndex(absId,iSupMod,iTower,iIphi,iIeta);
1016 //Gives SuperModule and Tower numbers
1017 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,
1018 iIphi, iIeta,iphi,ieta);
1019 //printf("Max id %d, iSM %d, col %d, row %d\n",absId,iSupMod,ieta,iphi);
1020 //printf("Max end---\n");
1023 //______________________________________
1024 void AliEMCALRecoUtils::InitParameters()
1026 // Initialize data members with default values
1028 fParticleType = kPhoton;
1029 fPosAlgo = kUnchanged;
1032 fNonLinearityFunction = kNoCorrection;
1033 fNonLinearThreshold = 30;
1035 fExoticCellFraction = 0.97;
1036 fExoticCellDiffTime = 1e6;
1037 fExoticCellMinAmplitude = 4.0;
1039 fAODFilterMask = 128;
1040 fAODHybridTracks = kFALSE;
1041 fAODTPCOnlyTracks = kTRUE;
1043 fCutEtaPhiSum = kTRUE;
1044 fCutEtaPhiSeparate = kFALSE;
1050 fClusterWindow = 100;
1055 fTrackCutsType = kLooseCut;
1058 fCutMinNClusterTPC = -1;
1059 fCutMinNClusterITS = -1;
1061 fCutMaxChi2PerClusterTPC = 1e10;
1062 fCutMaxChi2PerClusterITS = 1e10;
1064 fCutRequireTPCRefit = kFALSE;
1065 fCutRequireITSRefit = kFALSE;
1066 fCutAcceptKinkDaughters = kFALSE;
1068 fCutMaxDCAToVertexXY = 1e10;
1069 fCutMaxDCAToVertexZ = 1e10;
1070 fCutDCAToVertex2D = kFALSE;
1072 fCutRequireITSStandAlone = kFALSE; //MARCEL
1073 fCutRequireITSpureSA = kFALSE; //Marcel
1075 //Misalignment matrices
1076 for (Int_t i = 0; i < 15 ; i++)
1078 fMisalTransShift[i] = 0.;
1079 fMisalRotShift[i] = 0.;
1083 for (Int_t i = 0; i < 7 ; i++) fNonLinearityParams[i] = 0.;
1085 //For kBeamTestCorrectedv2 case, but default is no correction
1086 fNonLinearityParams[0] = 0.983504;
1087 fNonLinearityParams[1] = 0.210106;
1088 fNonLinearityParams[2] = 0.897274;
1089 fNonLinearityParams[3] = 0.0829064;
1090 fNonLinearityParams[4] = 152.299;
1091 fNonLinearityParams[5] = 31.5028;
1092 fNonLinearityParams[6] = 0.968;
1094 //Cluster energy smearing
1095 fSmearClusterEnergy = kFALSE;
1096 fSmearClusterParam[0] = 0.07; // * sqrt E term
1097 fSmearClusterParam[1] = 0.00; // * E term
1098 fSmearClusterParam[2] = 0.00; // constant
1101 //_____________________________________________________
1102 void AliEMCALRecoUtils::InitEMCALRecalibrationFactors()
1104 //Init EMCAL recalibration factors
1105 AliDebug(2,"AliCalorimeterUtils::InitEMCALRecalibrationFactors()");
1106 //In order to avoid rewriting the same histograms
1107 Bool_t oldStatus = TH1::AddDirectoryStatus();
1108 TH1::AddDirectory(kFALSE);
1110 fEMCALRecalibrationFactors = new TObjArray(12);
1111 for (int i = 0; i < 12; i++)
1112 fEMCALRecalibrationFactors->Add(new TH2F(Form("EMCALRecalFactors_SM%d",i),
1113 Form("EMCALRecalFactors_SM%d",i), 48, 0, 48, 24, 0, 24));
1114 //Init the histograms with 1
1115 for (Int_t sm = 0; sm < 12; sm++)
1117 for (Int_t i = 0; i < 48; i++)
1119 for (Int_t j = 0; j < 24; j++)
1121 SetEMCALChannelRecalibrationFactor(sm,i,j,1.);
1126 fEMCALRecalibrationFactors->SetOwner(kTRUE);
1127 fEMCALRecalibrationFactors->Compress();
1129 //In order to avoid rewriting the same histograms
1130 TH1::AddDirectory(oldStatus);
1133 //_________________________________________________________
1134 void AliEMCALRecoUtils::InitEMCALTimeRecalibrationFactors()
1136 //Init EMCAL recalibration factors
1137 AliDebug(2,"AliCalorimeterUtils::InitEMCALRecalibrationFactors()");
1138 //In order to avoid rewriting the same histograms
1139 Bool_t oldStatus = TH1::AddDirectoryStatus();
1140 TH1::AddDirectory(kFALSE);
1142 fEMCALTimeRecalibrationFactors = new TObjArray(4);
1143 for (int i = 0; i < 4; i++)
1144 fEMCALTimeRecalibrationFactors->Add(new TH1F(Form("hAllTimeAvBC%d",i),
1145 Form("hAllTimeAvBC%d",i),
1146 48*24*12,0.,48*24*12) );
1147 //Init the histograms with 1
1148 for (Int_t bc = 0; bc < 4; bc++)
1150 for (Int_t i = 0; i < 48*24*12; i++)
1151 SetEMCALChannelTimeRecalibrationFactor(bc,i,0.);
1154 fEMCALTimeRecalibrationFactors->SetOwner(kTRUE);
1155 fEMCALTimeRecalibrationFactors->Compress();
1157 //In order to avoid rewriting the same histograms
1158 TH1::AddDirectory(oldStatus);
1161 //____________________________________________________
1162 void AliEMCALRecoUtils::InitEMCALBadChannelStatusMap()
1164 //Init EMCAL bad channels map
1165 AliDebug(2,"AliEMCALRecoUtils::InitEMCALBadChannelStatusMap()");
1166 //In order to avoid rewriting the same histograms
1167 Bool_t oldStatus = TH1::AddDirectoryStatus();
1168 TH1::AddDirectory(kFALSE);
1170 fEMCALBadChannelMap = new TObjArray(12);
1171 //TH2F * hTemp = new TH2I("EMCALBadChannelMap","EMCAL SuperModule bad channel map", 48, 0, 48, 24, 0, 24);
1172 for (int i = 0; i < 12; i++)
1174 fEMCALBadChannelMap->Add(new TH2I(Form("EMCALBadChannelMap_Mod%d",i),Form("EMCALBadChannelMap_Mod%d",i), 48, 0, 48, 24, 0, 24));
1177 fEMCALBadChannelMap->SetOwner(kTRUE);
1178 fEMCALBadChannelMap->Compress();
1180 //In order to avoid rewriting the same histograms
1181 TH1::AddDirectory(oldStatus);
1184 //____________________________________________________________________________
1185 void AliEMCALRecoUtils::RecalibrateClusterEnergy(const AliEMCALGeometry* geom,
1186 AliVCluster * cluster,
1187 AliVCaloCells * cells,
1190 // Recalibrate the cluster energy and Time, considering the recalibration map
1191 // and the energy of the cells and time that compose the cluster.
1192 // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1195 AliInfo("Cluster pointer null!");
1199 //Get the cluster number of cells and list of absId, check what kind of cluster do we have.
1200 UShort_t * index = cluster->GetCellsAbsId() ;
1201 Double_t * fraction = cluster->GetCellsAmplitudeFraction() ;
1202 Int_t ncells = cluster->GetNCells();
1204 //Initialize some used variables
1207 Int_t icol =-1, irow =-1, imod=1;
1208 Float_t factor = 1, frac = 0;
1209 Int_t absIdMax = -1;
1212 //Loop on the cells, get the cell amplitude and recalibration factor, multiply and and to the new energy
1213 for (Int_t icell = 0; icell < ncells; icell++)
1215 absId = index[icell];
1216 frac = fraction[icell];
1217 if (frac < 1e-5) frac = 1; //in case of EMCAL, this is set as 0 since unfolding is off
1219 if (!fCellsRecalibrated && IsRecalibrationOn()) {
1221 Int_t iTower = -1, iIphi = -1, iIeta = -1;
1222 geom->GetCellIndex(absId,imod,iTower,iIphi,iIeta);
1223 if (fEMCALRecalibrationFactors->GetEntries() <= imod)
1225 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,irow,icol);
1226 factor = GetEMCALChannelRecalibrationFactor(imod,icol,irow);
1228 AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - recalibrate cell: module %d, col %d, row %d, cell fraction %f,recalibration factor %f, cell energy %f\n",
1229 imod,icol,irow,frac,factor,cells->GetCellAmplitude(absId)));
1233 energy += cells->GetCellAmplitude(absId)*factor*frac;
1235 if (emax < cells->GetCellAmplitude(absId)*factor*frac) {
1236 emax = cells->GetCellAmplitude(absId)*factor*frac;
1241 AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - Energy before %f, after %f \n",cluster->E(),energy));
1243 cluster->SetE(energy);
1245 // Recalculate time of cluster
1246 Double_t timeorg = cluster->GetTOF();
1248 Double_t time = cells->GetCellTime(absIdMax);
1249 if (!fCellsRecalibrated && IsTimeRecalibrationOn())
1250 RecalibrateCellTime(absIdMax,bc,time);
1252 cluster->SetTOF(time);
1254 AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - Time before %f, after %f \n",timeorg,cluster->GetTOF()));
1257 //_____________________________________________________________
1258 void AliEMCALRecoUtils::RecalibrateCells(AliVCaloCells * cells,
1261 // Recalibrate the cells time and energy, considering the recalibration map and the energy
1262 // of the cells that compose the cluster.
1263 // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1265 if (!IsRecalibrationOn() && !IsTimeRecalibrationOn() && !IsBadChannelsRemovalSwitchedOn())
1269 AliInfo("Cells pointer null!");
1274 Bool_t accept = kFALSE;
1277 Double_t ecellin = 0;
1278 Double_t tcellin = 0;
1282 Int_t nEMcell = cells->GetNumberOfCells() ;
1283 for (Int_t iCell = 0; iCell < nEMcell; iCell++)
1285 cells->GetCell( iCell, absId, ecellin, tcellin, mclabel, efrac );
1287 accept = AcceptCalibrateCell(absId, bc, ecell ,tcell ,cells);
1294 cells->SetCell(iCell,absId,ecell, tcell, mclabel, efrac);
1297 fCellsRecalibrated = kTRUE;
1300 //_______________________________________________________________________________________________________
1301 void AliEMCALRecoUtils::RecalibrateCellTime(Int_t absId, Int_t bc, Double_t & celltime) const
1303 // Recalibrate time of cell with absID considering the recalibration map
1304 // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1306 if (!fCellsRecalibrated && IsTimeRecalibrationOn() && bc >= 0) {
1307 celltime -= GetEMCALChannelTimeRecalibrationFactor(bc%4,absId)*1.e-9; ;
1311 //______________________________________________________________________________
1312 void AliEMCALRecoUtils::RecalculateClusterPosition(const AliEMCALGeometry *geom,
1313 AliVCaloCells* cells,
1316 //For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1319 AliInfo("Cluster pointer null!");
1323 if (fPosAlgo==kPosTowerGlobal) RecalculateClusterPositionFromTowerGlobal( geom, cells, clu);
1324 else if (fPosAlgo==kPosTowerIndex) RecalculateClusterPositionFromTowerIndex ( geom, cells, clu);
1325 else AliDebug(2,"Algorithm to recalculate position not selected, do nothing.");
1328 //_____________________________________________________________________________________________
1329 void AliEMCALRecoUtils::RecalculateClusterPositionFromTowerGlobal(const AliEMCALGeometry *geom,
1330 AliVCaloCells* cells,
1333 // For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1334 // The algorithm is a copy of what is done in AliEMCALRecPoint
1336 Double_t eCell = 0.;
1337 Float_t fraction = 1.;
1338 Float_t recalFactor = 1.;
1341 Int_t iTower = -1, iIphi = -1, iIeta = -1;
1342 Int_t iSupModMax = -1, iSM=-1, iphi = -1, ieta = -1;
1343 Float_t weight = 0., totalWeight=0.;
1344 Float_t newPos[3] = {0,0,0};
1345 Double_t pLocal[3], pGlobal[3];
1346 Bool_t shared = kFALSE;
1348 Float_t clEnergy = clu->E(); //Energy already recalibrated previously
1351 GetMaxEnergyCell(geom, cells, clu, absId, iSupModMax, ieta, iphi,shared);
1352 Double_t depth = GetDepth(clEnergy,fParticleType,iSupModMax) ;
1354 //printf("** Cluster energy %f, ncells %d, depth %f\n",clEnergy,clu->GetNCells(),depth);
1356 for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++)
1358 absId = clu->GetCellAbsId(iDig);
1359 fraction = clu->GetCellAmplitudeFraction(iDig);
1360 if (fraction < 1e-4) fraction = 1.; // in case unfolding is off
1362 if (!fCellsRecalibrated) {
1363 geom->GetCellIndex(absId,iSM,iTower,iIphi,iIeta);
1364 geom->GetCellPhiEtaIndexInSModule(iSM,iTower,iIphi, iIeta,iphi,ieta);
1365 if (IsRecalibrationOn()) {
1366 recalFactor = GetEMCALChannelRecalibrationFactor(iSM,ieta,iphi);
1370 eCell = cells->GetCellAmplitude(absId)*fraction*recalFactor;
1372 weight = GetCellWeight(eCell,clEnergy);
1373 totalWeight += weight;
1375 geom->RelPosCellInSModule(absId,depth,pLocal[0],pLocal[1],pLocal[2]);
1376 //printf("pLocal (%f,%f,%f), SM %d, absId %d\n",pLocal[0],pLocal[1],pLocal[2],iSupModMax,absId);
1377 geom->GetGlobal(pLocal,pGlobal,iSupModMax);
1378 //printf("pLocal (%f,%f,%f)\n",pGlobal[0],pGlobal[1],pGlobal[2]);
1380 for (int i=0; i<3; i++ ) newPos[i] += (weight*pGlobal[i]);
1383 if (totalWeight>0) {
1384 for (int i=0; i<3; i++ ) newPos[i] /= totalWeight;
1387 //Float_t pos[]={0,0,0};
1388 //clu->GetPosition(pos);
1389 //printf("OldPos : %2.3f,%2.3f,%2.3f\n",pos[0],pos[1],pos[2]);
1390 //printf("NewPos : %2.3f,%2.3f,%2.3f\n",newPos[0],newPos[1],newPos[2]);
1392 if (iSupModMax > 1) { //sector 1
1393 newPos[0] +=fMisalTransShift[3];//-=3.093;
1394 newPos[1] +=fMisalTransShift[4];//+=6.82;
1395 newPos[2] +=fMisalTransShift[5];//+=1.635;
1396 //printf(" + : %2.3f,%2.3f,%2.3f\n",fMisalTransShift[3],fMisalTransShift[4],fMisalTransShift[5]);
1398 newPos[0] +=fMisalTransShift[0];//+=1.134;
1399 newPos[1] +=fMisalTransShift[1];//+=8.2;
1400 newPos[2] +=fMisalTransShift[2];//+=1.197;
1401 //printf(" + : %2.3f,%2.3f,%2.3f\n",fMisalTransShift[0],fMisalTransShift[1],fMisalTransShift[2]);
1403 //printf("NewPos : %2.3f,%2.3f,%2.3f\n",newPos[0],newPos[1],newPos[2]);
1405 clu->SetPosition(newPos);
1408 //____________________________________________________________________________________________
1409 void AliEMCALRecoUtils::RecalculateClusterPositionFromTowerIndex(const AliEMCALGeometry *geom,
1410 AliVCaloCells* cells,
1413 // For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1414 // The algorithm works with the tower indeces, averages the indeces and from them it calculates the global position
1416 Double_t eCell = 1.;
1417 Float_t fraction = 1.;
1418 Float_t recalFactor = 1.;
1422 Int_t iIphi = -1, iIeta = -1;
1423 Int_t iSupMod = -1, iSupModMax = -1;
1424 Int_t iphi = -1, ieta =-1;
1425 Bool_t shared = kFALSE;
1427 Float_t clEnergy = clu->E(); //Energy already recalibrated previously.
1431 GetMaxEnergyCell(geom, cells, clu, absId, iSupModMax, ieta, iphi,shared);
1432 Float_t depth = GetDepth(clEnergy,fParticleType,iSupMod) ;
1434 Float_t weight = 0., weightedCol = 0., weightedRow = 0., totalWeight=0.;
1435 Bool_t areInSameSM = kTRUE; //exclude clusters with cells in different SMs for now
1436 Int_t startingSM = -1;
1438 for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++)
1440 absId = clu->GetCellAbsId(iDig);
1441 fraction = clu->GetCellAmplitudeFraction(iDig);
1442 if (fraction < 1e-4) fraction = 1.; // in case unfolding is off
1444 if (iDig==0) startingSM = iSupMod;
1445 else if (iSupMod != startingSM) areInSameSM = kFALSE;
1447 eCell = cells->GetCellAmplitude(absId);
1449 geom->GetCellIndex(absId,iSupMod,iTower,iIphi,iIeta);
1450 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);
1452 if (!fCellsRecalibrated)
1454 if (IsRecalibrationOn()) {
1455 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1459 eCell = cells->GetCellAmplitude(absId)*fraction*recalFactor;
1461 weight = GetCellWeight(eCell,clEnergy);
1462 if (weight < 0) weight = 0;
1463 totalWeight += weight;
1464 weightedCol += ieta*weight;
1465 weightedRow += iphi*weight;
1467 //printf("Max cell? cell %d, amplitude org %f, fraction %f, recalibration %f, amplitude new %f \n",cellAbsId, cells->GetCellAmplitude(cellAbsId), fraction, recalFactor, eCell) ;
1470 Float_t xyzNew[]={0.,0.,0.};
1471 if (areInSameSM == kTRUE) {
1472 //printf("In Same SM\n");
1473 weightedCol = weightedCol/totalWeight;
1474 weightedRow = weightedRow/totalWeight;
1475 geom->RecalculateTowerPosition(weightedRow, weightedCol, iSupModMax, depth, fMisalTransShift, fMisalRotShift, xyzNew);
1479 //printf("In Different SM\n");
1480 geom->RecalculateTowerPosition(iphi, ieta, iSupModMax, depth, fMisalTransShift, fMisalRotShift, xyzNew);
1483 clu->SetPosition(xyzNew);
1486 //___________________________________________________________________________________________
1487 void AliEMCALRecoUtils::RecalculateClusterDistanceToBadChannel(const AliEMCALGeometry * geom,
1488 AliVCaloCells* cells,
1489 AliVCluster * cluster)
1491 //re-evaluate distance to bad channel with updated bad map
1493 if (!fRecalDistToBadChannels) return;
1497 AliInfo("Cluster pointer null!");
1501 //Get channels map of the supermodule where the cluster is.
1502 Int_t absIdMax = -1, iSupMod =-1, icolM = -1, irowM = -1;
1503 Bool_t shared = kFALSE;
1504 GetMaxEnergyCell(geom, cells, cluster, absIdMax, iSupMod, icolM, irowM, shared);
1505 TH2D* hMap = (TH2D*)fEMCALBadChannelMap->At(iSupMod);
1508 Float_t minDist = 10000.;
1511 //Loop on tower status map
1512 for (Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++)
1514 for (Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++)
1516 //Check if tower is bad.
1517 if (hMap->GetBinContent(icol,irow)==0) continue;
1518 //printf("AliEMCALRecoUtils::RecalculateDistanceToBadChannels() - \n \t Bad channel in SM %d, col %d, row %d, \n \t Cluster max in col %d, row %d\n",
1519 // iSupMod,icol, irow, icolM,irowM);
1521 dRrow=TMath::Abs(irowM-irow);
1522 dRcol=TMath::Abs(icolM-icol);
1523 dist=TMath::Sqrt(dRrow*dRrow+dRcol*dRcol);
1526 //printf("MIN DISTANCE TO BAD %2.2f\n",dist);
1532 //In case the cluster is shared by 2 SuperModules, need to check the map of the second Super Module
1536 Int_t iSupMod2 = -1;
1538 //The only possible combinations are (0,1), (2,3) ... (8,9)
1539 if (iSupMod%2) iSupMod2 = iSupMod-1;
1540 else iSupMod2 = iSupMod+1;
1541 hMap2 = (TH2D*)fEMCALBadChannelMap->At(iSupMod2);
1543 //Loop on tower status map of second super module
1544 for (Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++)
1546 for (Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++)
1548 //Check if tower is bad.
1549 if (hMap2->GetBinContent(icol,irow)==0)
1551 //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",
1552 // iSupMod2,icol, irow,iSupMod,icolM,irowM);
1553 dRrow=TMath::Abs(irow-irowM);
1556 dRcol=TMath::Abs(icol-(AliEMCALGeoParams::fgkEMCALCols+icolM));
1558 dRcol=TMath::Abs(AliEMCALGeoParams::fgkEMCALCols+icol-icolM);
1561 dist=TMath::Sqrt(dRrow*dRrow+dRcol*dRcol);
1562 if (dist < minDist) minDist = dist;
1565 }// shared cluster in 2 SuperModules
1567 AliDebug(2,Form("Max cluster cell (SM,col,row)=(%d %d %d) - Distance to Bad Channel %2.2f",iSupMod, icolM, irowM, minDist));
1568 cluster->SetDistanceToBadChannel(minDist);
1571 //__________________________________________________________________
1572 void AliEMCALRecoUtils::RecalculateClusterPID(AliVCluster * cluster)
1574 //re-evaluate identification parameters with bayesian
1577 AliInfo("Cluster pointer null!");
1581 if (cluster->GetM02() != 0)
1582 fPIDUtils->ComputePID(cluster->E(),cluster->GetM02());
1584 Float_t pidlist[AliPID::kSPECIESCN+1];
1585 for (Int_t i = 0; i < AliPID::kSPECIESCN+1; i++) pidlist[i] = fPIDUtils->GetPIDFinal(i);
1587 cluster->SetPID(pidlist);
1590 //___________________________________________________________________________________________________________________
1591 void AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(const AliEMCALGeometry * geom,
1592 AliVCaloCells* cells,
1593 AliVCluster * cluster,
1594 Float_t & l0, Float_t & l1,
1595 Float_t & disp, Float_t & dEta, Float_t & dPhi,
1596 Float_t & sEta, Float_t & sPhi, Float_t & sEtaPhi)
1598 // Calculates new center of gravity in the local EMCAL-module coordinates
1599 // and tranfers into global ALICE coordinates
1600 // Calculates Dispersion and main axis
1603 AliInfo("Cluster pointer null!");
1607 Double_t eCell = 0.;
1608 Float_t fraction = 1.;
1609 Float_t recalFactor = 1.;
1617 Double_t etai = -1.;
1618 Double_t phii = -1.;
1623 Double_t etaMean = 0.;
1624 Double_t phiMean = 0.;
1626 //Loop on cells, calculate the cluster energy, in case a cut on cell energy is added
1627 // and to check if the cluster is between 2 SM in eta
1629 Bool_t shared = kFALSE;
1632 for (Int_t iDigit=0; iDigit < cluster->GetNCells(); iDigit++)
1634 //Get from the absid the supermodule, tower and eta/phi numbers
1635 geom->GetCellIndex(cluster->GetCellAbsId(iDigit),iSupMod,iTower,iIphi,iIeta);
1636 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi,iIeta, iphi,ieta);
1638 //Check if there are cells of different SM
1639 if (iDigit == 0 ) iSM0 = iSupMod;
1640 else if (iSupMod!= iSM0) shared = kTRUE;
1642 //Get the cell energy, if recalibration is on, apply factors
1643 fraction = cluster->GetCellAmplitudeFraction(iDigit);
1644 if (fraction < 1e-4) fraction = 1.; // in case unfolding is off
1646 if (IsRecalibrationOn()) {
1647 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1650 eCell = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1657 for (Int_t iDigit=0; iDigit < cluster->GetNCells(); iDigit++)
1659 //Get from the absid the supermodule, tower and eta/phi numbers
1660 geom->GetCellIndex(cluster->GetCellAbsId(iDigit),iSupMod,iTower,iIphi,iIeta);
1661 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi,iIeta, iphi,ieta);
1663 //Get the cell energy, if recalibration is on, apply factors
1664 fraction = cluster->GetCellAmplitudeFraction(iDigit);
1665 if (fraction < 1e-4) fraction = 1.; // in case unfolding is off
1667 if (!fCellsRecalibrated) {
1668 if (IsRecalibrationOn()) {
1669 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1673 eCell = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1675 // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
1676 // C Side impair SM, nSupMod%2=1; A side pair SM, nSupMod%2=0
1677 if (shared && iSupMod%2) ieta+=AliEMCALGeoParams::fgkEMCALCols;
1679 if (cluster->E() > 0 && eCell > 0) {
1680 w = GetCellWeight(eCell,cluster->E());
1682 etai=(Double_t)ieta;
1683 phii=(Double_t)iphi;
1689 sEta += w * etai * etai ;
1690 etaMean += w * etai ;
1691 sPhi += w * phii * phii ;
1692 phiMean += w * phii ;
1693 sEtaPhi += w * etai * phii ;
1696 AliError(Form("Wrong energy %f and/or amplitude %f\n", eCell, cluster->E()));
1699 //Normalize to the weight
1704 AliError(Form("Wrong weight %f\n", wtot));
1706 //Calculate dispersion
1707 for (Int_t iDigit=0; iDigit < cluster->GetNCells(); iDigit++)
1709 //Get from the absid the supermodule, tower and eta/phi numbers
1710 geom->GetCellIndex(cluster->GetCellAbsId(iDigit),iSupMod,iTower,iIphi,iIeta);
1711 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi,iIeta, iphi,ieta);
1713 //Get the cell energy, if recalibration is on, apply factors
1714 fraction = cluster->GetCellAmplitudeFraction(iDigit);
1715 if (fraction < 1e-4) fraction = 1.; // in case unfolding is off
1716 if (IsRecalibrationOn()) {
1717 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1719 eCell = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1721 // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
1722 // C Side impair SM, nSupMod%2=1; A side pair SM, nSupMod%2=0
1723 if (shared && iSupMod%2) ieta+=AliEMCALGeoParams::fgkEMCALCols;
1725 if (cluster->E() > 0 && eCell > 0) {
1726 w = GetCellWeight(eCell,cluster->E());
1728 etai=(Double_t)ieta;
1729 phii=(Double_t)iphi;
1731 disp += w *((etai-etaMean)*(etai-etaMean)+(phii-phiMean)*(phii-phiMean));
1732 dEta += w * (etai-etaMean)*(etai-etaMean) ;
1733 dPhi += w * (phii-phiMean)*(phii-phiMean) ;
1736 AliError(Form("Wrong energy %f and/or amplitude %f\n", eCell, cluster->E()));
1739 //Normalize to the weigth and set shower shape parameters
1740 if (wtot > 0 && nstat > 1) {
1748 sEta -= etaMean * etaMean ;
1749 sPhi -= phiMean * phiMean ;
1750 sEtaPhi -= etaMean * phiMean ;
1752 l0 = (0.5 * (sEta + sPhi) + TMath::Sqrt( 0.25 * (sEta - sPhi) * (sEta - sPhi) + sEtaPhi * sEtaPhi ));
1753 l1 = (0.5 * (sEta + sPhi) - TMath::Sqrt( 0.25 * (sEta - sPhi) * (sEta - sPhi) + sEtaPhi * sEtaPhi ));
1757 dEta = 0. ; dPhi = 0. ; disp = 0. ;
1758 sEta = 0. ; sPhi = 0. ; sEtaPhi = 0. ;
1762 //____________________________________________________________________________________________
1763 void AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(const AliEMCALGeometry * geom,
1764 AliVCaloCells* cells,
1765 AliVCluster * cluster)
1767 // Calculates new center of gravity in the local EMCAL-module coordinates
1768 // and tranfers into global ALICE coordinates
1769 // Calculates Dispersion and main axis and puts them into the cluster
1771 Float_t l0 = 0., l1 = 0.;
1772 Float_t disp = 0., dEta = 0., dPhi = 0.;
1773 Float_t sEta = 0., sPhi = 0., sEtaPhi = 0.;
1775 AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(geom,cells,cluster,l0,l1,disp,
1776 dEta, dPhi, sEta, sPhi, sEtaPhi);
1778 cluster->SetM02(l0);
1779 cluster->SetM20(l1);
1780 if (disp > 0. ) cluster->SetDispersion(TMath::Sqrt(disp)) ;
1784 //____________________________________________________________________________
1785 void AliEMCALRecoUtils::FindMatches(AliVEvent *event,
1786 TObjArray * clusterArr,
1787 const AliEMCALGeometry *geom)
1789 //This function should be called before the cluster loop
1790 //Before call this function, please recalculate the cluster positions
1791 //Given the input event, loop over all the tracks, select the closest cluster as matched with fCutR
1792 //Store matched cluster indexes and residuals
1794 fMatchedTrackIndex ->Reset();
1795 fMatchedClusterIndex->Reset();
1796 fResidualPhi->Reset();
1797 fResidualEta->Reset();
1799 fMatchedTrackIndex ->Set(1000);
1800 fMatchedClusterIndex->Set(1000);
1801 fResidualPhi->Set(1000);
1802 fResidualEta->Set(1000);
1804 AliESDEvent* esdevent = dynamic_cast<AliESDEvent*> (event);
1805 AliAODEvent* aodevent = dynamic_cast<AliAODEvent*> (event);
1807 // init the magnetic field if not already on
1808 if (!TGeoGlobalMagField::Instance()->GetField()) {
1809 if (!event->InitMagneticField())
1811 AliInfo("Mag Field not initialized, null esd/aod evetn pointers");
1816 UInt_t mask1 = esdevent->GetESDRun()->GetDetectorsInDAQ();
1817 UInt_t mask2 = esdevent->GetESDRun()->GetDetectorsInReco();
1818 Bool_t desc1 = (mask1 >> 3) & 0x1;
1819 Bool_t desc2 = (mask2 >> 3) & 0x1;
1820 if (desc1==0 || desc2==0) {
1821 // AliError(Form("TPC not in DAQ/RECO: %u (%u)/%u (%u)",
1822 // mask1, esdevent->GetESDRun()->GetDetectorsInReco(),
1823 // mask2, esdevent->GetESDRun()->GetDetectorsInDAQ()));
1828 TObjArray *clusterArray = 0x0;
1830 clusterArray = new TObjArray(event->GetNumberOfCaloClusters());
1831 for (Int_t icl=0; icl<event->GetNumberOfCaloClusters(); icl++)
1833 AliVCluster *cluster = (AliVCluster*) event->GetCaloCluster(icl);
1834 if (geom && !IsGoodCluster(cluster,geom,(AliVCaloCells*)event->GetEMCALCells()))
1836 clusterArray->AddAt(cluster,icl);
1842 for (Int_t i=0; i<21;i++) cv[i]=0;
1843 for (Int_t itr=0; itr<event->GetNumberOfTracks(); itr++)
1845 AliExternalTrackParam *trackParam = 0;
1847 //If the input event is ESD, the starting point for extrapolation is TPCOut, if available, or TPCInner
1848 AliESDtrack *esdTrack = 0;
1849 AliAODTrack *aodTrack = 0;
1851 esdTrack = esdevent->GetTrack(itr);
1852 if (!esdTrack) continue;
1853 if (!IsAccepted(esdTrack)) continue;
1854 if (esdTrack->Pt()<fCutMinTrackPt) continue;
1855 Double_t phi = esdTrack->Phi()*TMath::RadToDeg();
1856 if (TMath::Abs(esdTrack->Eta())>0.9 || phi <= 10 || phi >= 250 ) continue;
1858 trackParam = const_cast<AliExternalTrackParam*>(esdTrack->GetInnerParam()); // if TPC Available
1860 trackParam = new AliExternalTrackParam(*esdTrack); // If ITS Track Standing alone
1863 //If the input event is AOD, the starting point for extrapolation is at vertex
1864 //AOD tracks are selected according to its filterbit.
1865 else if (aodevent) {
1866 aodTrack = dynamic_cast<AliAODTrack*>(aodevent->GetTrack(itr));
1867 if(!aodTrack) AliFatal("Not a standard AOD");
1868 if (!aodTrack) continue;
1870 if (fAODTPCOnlyTracks) { // Match with TPC only tracks, default from May 2013, before filter bit 32
1871 //printf("Match with TPC only tracks, accept? %d, test bit 128 <%d> \n", aodTrack->IsTPCOnly(), aodTrack->TestFilterMask(128));
1872 if (!aodTrack->IsTPCConstrained()) continue ;
1873 } else if (fAODHybridTracks) { // Match with hybrid tracks
1874 //printf("Match with Hybrid tracks, accept? %d \n", aodTrack->IsHybridGlobalConstrainedGlobal());
1875 if (!aodTrack->IsHybridGlobalConstrainedGlobal()) continue ;
1876 } else { // Match with tracks on a mask
1877 //printf("Match with tracks having filter bit mask %d, accept? %d \n",fAODFilterMask,aodTrack->TestFilterMask(fAODFilterMask));
1878 if (!aodTrack->TestFilterMask(fAODFilterMask) ) continue; //Select AOD tracks
1881 if (aodTrack->Pt()<fCutMinTrackPt) continue;
1883 Double_t phi = aodTrack->Phi()*TMath::RadToDeg();
1884 if (TMath::Abs(aodTrack->Eta())>0.9 || phi <= 10 || phi >= 250 )
1886 Double_t pos[3],mom[3];
1887 aodTrack->GetXYZ(pos);
1888 aodTrack->GetPxPyPz(mom);
1889 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()));
1891 trackParam= new AliExternalTrackParam(pos,mom,cv,aodTrack->Charge());
1894 //Return if the input data is not "AOD" or "ESD"
1896 printf("Wrong input data type! Should be \"AOD\" or \"ESD\"\n");
1898 clusterArray->Clear();
1899 delete clusterArray;
1904 if (!trackParam) continue;
1906 //Extrapolate the track to EMCal surface
1907 AliExternalTrackParam emcalParam(*trackParam);
1908 Float_t eta, phi, pt;
1909 if (!ExtrapolateTrackToEMCalSurface(&emcalParam, fEMCalSurfaceDistance, fMass, fStepSurface, eta, phi, pt)) {
1910 if (aodevent && trackParam) delete trackParam;
1911 if (fITSTrackSA && trackParam) delete trackParam;
1915 if (TMath::Abs(eta)>0.75 || (phi) < 70*TMath::DegToRad() || (phi) > 190*TMath::DegToRad()) {
1916 if (aodevent && trackParam) delete trackParam;
1917 if (fITSTrackSA && trackParam) delete trackParam;
1921 //Find matched clusters
1923 Float_t dEta = -999, dPhi = -999;
1925 index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArray, dEta, dPhi);
1927 index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArr, dEta, dPhi);
1931 fMatchedTrackIndex ->AddAt(itr,matched);
1932 fMatchedClusterIndex ->AddAt(index,matched);
1933 fResidualEta ->AddAt(dEta,matched);
1934 fResidualPhi ->AddAt(dPhi,matched);
1937 if (aodevent && trackParam) delete trackParam;
1938 if (fITSTrackSA && trackParam) delete trackParam;
1942 clusterArray->Clear();
1943 delete clusterArray;
1946 AliDebug(2,Form("Number of matched pairs = %d !\n",matched));
1948 fMatchedTrackIndex ->Set(matched);
1949 fMatchedClusterIndex ->Set(matched);
1950 fResidualPhi ->Set(matched);
1951 fResidualEta ->Set(matched);
1954 //________________________________________________________________________________
1955 Int_t AliEMCALRecoUtils::FindMatchedClusterInEvent(const AliESDtrack *track,
1956 const AliVEvent *event,
1957 const AliEMCALGeometry *geom,
1958 Float_t &dEta, Float_t &dPhi)
1961 // This function returns the index of matched cluster to input track
1962 // Returns -1 if no match is found
1964 Double_t phiV = track->Phi()*TMath::RadToDeg();
1965 if (TMath::Abs(track->Eta())>0.9 || phiV <= 10 || phiV >= 250 ) return index;
1966 AliExternalTrackParam *trackParam = 0;
1968 trackParam = const_cast<AliExternalTrackParam*>(track->GetInnerParam()); // If TPC
1970 trackParam = new AliExternalTrackParam(*track);
1972 if (!trackParam) return index;
1973 AliExternalTrackParam emcalParam(*trackParam);
1974 Float_t eta, phi, pt;
1976 if (!ExtrapolateTrackToEMCalSurface(&emcalParam, fEMCalSurfaceDistance, fMass, fStepSurface, eta, phi, pt)) {
1977 if (fITSTrackSA) delete trackParam;
1980 if (TMath::Abs(eta)>0.75 || (phi) < 70*TMath::DegToRad() || (phi) > 190*TMath::DegToRad()) {
1981 if (fITSTrackSA) delete trackParam;
1985 TObjArray *clusterArr = new TObjArray(event->GetNumberOfCaloClusters());
1987 for (Int_t icl=0; icl<event->GetNumberOfCaloClusters(); icl++)
1989 AliVCluster *cluster = (AliVCluster*) event->GetCaloCluster(icl);
1990 if (geom && !IsGoodCluster(cluster,geom,(AliVCaloCells*)event->GetEMCALCells())) continue;
1991 clusterArr->AddAt(cluster,icl);
1994 index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArr, dEta, dPhi);
1995 clusterArr->Clear();
1997 if (fITSTrackSA) delete trackParam;
2002 //_______________________________________________________________________________________________
2003 Int_t AliEMCALRecoUtils::FindMatchedClusterInClusterArr(const AliExternalTrackParam *emcalParam,
2004 AliExternalTrackParam *trkParam,
2005 const TObjArray * clusterArr,
2006 Float_t &dEta, Float_t &dPhi)
2008 // Find matched cluster in array
2010 dEta=-999, dPhi=-999;
2011 Float_t dRMax = fCutR, dEtaMax=fCutEta, dPhiMax=fCutPhi;
2013 Float_t tmpEta=-999, tmpPhi=-999;
2015 Double_t exPos[3] = {0.,0.,0.};
2016 if (!emcalParam->GetXYZ(exPos)) return index;
2018 Float_t clsPos[3] = {0.,0.,0.};
2019 for (Int_t icl=0; icl<clusterArr->GetEntriesFast(); icl++)
2021 AliVCluster *cluster = dynamic_cast<AliVCluster*> (clusterArr->At(icl)) ;
2022 if (!cluster || !cluster->IsEMCAL()) continue;
2023 cluster->GetPosition(clsPos);
2024 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));
2025 if (dR > fClusterWindow) continue;
2027 AliExternalTrackParam trkPamTmp (*trkParam);//Retrieve the starting point every time before the extrapolation
2028 if (!ExtrapolateTrackToCluster(&trkPamTmp, cluster, fMass, fStepCluster, tmpEta, tmpPhi)) continue;
2029 if (fCutEtaPhiSum) {
2030 Float_t tmpR=TMath::Sqrt(tmpEta*tmpEta + tmpPhi*tmpPhi);
2037 } else if (fCutEtaPhiSeparate) {
2038 if (TMath::Abs(tmpEta)<TMath::Abs(dEtaMax) && TMath::Abs(tmpPhi)<TMath::Abs(dPhiMax)) {
2044 printf("Error: please specify your cut criteria\n");
2045 printf("To cut on sqrt(dEta^2+dPhi^2), use: SwitchOnCutEtaPhiSum()\n");
2046 printf("To cut on dEta and dPhi separately, use: SwitchOnCutEtaPhiSeparate()\n");
2057 //------------------------------------------------------------------------------------
2058 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToEMCalSurface(AliVTrack *track,
2059 Double_t emcalR, Double_t mass,
2060 Double_t step, Double_t minpt,
2061 Bool_t useMassForTracking)
2063 // Extrapolate track to EMCAL surface
2065 track->SetTrackPhiEtaPtOnEMCal(-999, -999, -999);
2067 if (track->Pt()<minpt)
2070 Double_t phi = track->Phi()*TMath::RadToDeg();
2071 if (TMath::Abs(track->Eta())>0.9 || phi <= 10 || phi >= 250)
2074 AliESDtrack *esdt = dynamic_cast<AliESDtrack*>(track);
2075 AliAODTrack *aodt = 0;
2077 aodt = dynamic_cast<AliAODTrack*>(track);
2082 // Select the mass hypothesis
2085 Bool_t onlyTPC = kFALSE;
2086 if ( mass == -99 ) onlyTPC=kTRUE;
2090 if ( useMassForTracking ) mass = esdt->GetMassForTracking();
2091 else mass = esdt->GetMass(onlyTPC);
2095 if ( useMassForTracking ) mass = aodt->GetMassForTracking();
2096 else mass = aodt->M();
2100 AliExternalTrackParam *trackParam = 0;
2102 const AliExternalTrackParam *in = esdt->GetInnerParam();
2105 trackParam = new AliExternalTrackParam(*in);
2107 Double_t xyz[3] = {0}, pxpypz[3] = {0}, cv[21] = {0};
2108 aodt->PxPyPz(pxpypz);
2110 aodt->GetCovarianceXYZPxPyPz(cv);
2111 trackParam = new AliExternalTrackParam(xyz,pxpypz,cv,aodt->Charge());
2116 Float_t etaout=-999, phiout=-999, ptout=-999;
2117 Bool_t ret = ExtrapolateTrackToEMCalSurface(trackParam,
2127 if (TMath::Abs(etaout)>0.75 || (phiout<70*TMath::DegToRad()) || (phiout>190*TMath::DegToRad()))
2129 track->SetTrackPhiEtaPtOnEMCal(phiout, etaout, ptout);
2134 //------------------------------------------------------------------------------------
2135 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToEMCalSurface(AliExternalTrackParam *trkParam,
2143 //Extrapolate track to EMCAL surface
2145 eta = -999, phi = -999, pt = -999;
2146 if (!trkParam) return kFALSE;
2147 if (!AliTrackerBase::PropagateTrackToBxByBz(trkParam, emcalR, mass, step, kTRUE, 0.8, -1)) return kFALSE;
2148 Double_t trkPos[3] = {0.,0.,0.};
2149 if (!trkParam->GetXYZ(trkPos)) return kFALSE;
2150 TVector3 trkPosVec(trkPos[0],trkPos[1],trkPos[2]);
2151 eta = trkPosVec.Eta();
2152 phi = trkPosVec.Phi();
2153 pt = trkParam->Pt();
2155 phi += 2*TMath::Pi();
2160 //-----------------------------------------------------------------------------------
2161 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToPosition(AliExternalTrackParam *trkParam,
2162 const Float_t *clsPos,
2169 //Return the residual by extrapolating a track param to a global position
2173 if (!trkParam) return kFALSE;
2174 Double_t trkPos[3] = {0.,0.,0.};
2175 TVector3 vec(clsPos[0],clsPos[1],clsPos[2]);
2176 Double_t alpha = ((int)(vec.Phi()*TMath::RadToDeg()/20)+0.5)*20*TMath::DegToRad();
2177 vec.RotateZ(-alpha); //Rotate the cluster to the local extrapolation coordinate system
2178 if (!AliTrackerBase::PropagateTrackToBxByBz(trkParam, vec.X(), mass, step,kTRUE, 0.8, -1)) return kFALSE;
2179 if (!trkParam->GetXYZ(trkPos)) return kFALSE; //Get the extrapolated global position
2181 TVector3 clsPosVec(clsPos[0],clsPos[1],clsPos[2]);
2182 TVector3 trkPosVec(trkPos[0],trkPos[1],trkPos[2]);
2184 // track cluster matching
2185 tmpPhi = clsPosVec.DeltaPhi(trkPosVec); // tmpPhi is between -pi and pi
2186 tmpEta = clsPosVec.Eta()-trkPosVec.Eta();
2191 //----------------------------------------------------------------------------------
2192 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToCluster(AliExternalTrackParam *trkParam,
2193 const AliVCluster *cluster,
2200 //Return the residual by extrapolating a track param to a cluster
2204 if (!cluster || !trkParam)
2207 Float_t clsPos[3] = {0.,0.,0.};
2208 cluster->GetPosition(clsPos);
2210 return ExtrapolateTrackToPosition(trkParam, clsPos, mass, step, tmpEta, tmpPhi);
2213 //---------------------------------------------------------------------------------
2214 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToCluster(AliExternalTrackParam *trkParam,
2215 const AliVCluster *cluster,
2220 //Return the residual by extrapolating a track param to a clusterfStepCluster
2223 return ExtrapolateTrackToCluster(trkParam, cluster, fMass, fStepCluster, tmpEta, tmpPhi);
2226 //_______________________________________________________________________
2227 void AliEMCALRecoUtils::GetMatchedResiduals(Int_t clsIndex,
2228 Float_t &dEta, Float_t &dPhi)
2230 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2231 //Get the residuals dEta and dPhi for this cluster to the closest track
2232 //Works with ESDs and AODs
2234 if (FindMatchedPosForCluster(clsIndex) >= 999) {
2235 AliDebug(2,"No matched tracks found!\n");
2240 dEta = fResidualEta->At(FindMatchedPosForCluster(clsIndex));
2241 dPhi = fResidualPhi->At(FindMatchedPosForCluster(clsIndex));
2244 //______________________________________________________________________________________________
2245 void AliEMCALRecoUtils::GetMatchedClusterResiduals(Int_t trkIndex, Float_t &dEta, Float_t &dPhi)
2247 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2248 //Get the residuals dEta and dPhi for this track to the closest cluster
2249 //Works with ESDs and AODs
2251 if (FindMatchedPosForTrack(trkIndex) >= 999) {
2252 AliDebug(2,"No matched cluster found!\n");
2257 dEta = fResidualEta->At(FindMatchedPosForTrack(trkIndex));
2258 dPhi = fResidualPhi->At(FindMatchedPosForTrack(trkIndex));
2261 //__________________________________________________________
2262 Int_t AliEMCALRecoUtils::GetMatchedTrackIndex(Int_t clsIndex)
2264 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2265 //Get the index of matched track to this cluster
2266 //Works with ESDs and AODs
2268 if (IsClusterMatched(clsIndex))
2269 return fMatchedTrackIndex->At(FindMatchedPosForCluster(clsIndex));
2274 //__________________________________________________________
2275 Int_t AliEMCALRecoUtils::GetMatchedClusterIndex(Int_t trkIndex)
2277 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2278 //Get the index of matched cluster to this track
2279 //Works with ESDs and AODs
2281 if (IsTrackMatched(trkIndex))
2282 return fMatchedClusterIndex->At(FindMatchedPosForTrack(trkIndex));
2287 //______________________________________________________________
2288 Bool_t AliEMCALRecoUtils::IsClusterMatched(Int_t clsIndex) const
2290 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2291 //Returns if the cluster has a match
2292 if (FindMatchedPosForCluster(clsIndex) < 999)
2298 //____________________________________________________________
2299 Bool_t AliEMCALRecoUtils::IsTrackMatched(Int_t trkIndex) const
2301 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2302 //Returns if the track has a match
2303 if (FindMatchedPosForTrack(trkIndex) < 999)
2309 //______________________________________________________________________
2310 UInt_t AliEMCALRecoUtils::FindMatchedPosForCluster(Int_t clsIndex) const
2312 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2313 //Returns the position of the match in the fMatchedClusterIndex array
2314 Float_t tmpR = fCutR;
2317 for (Int_t i=0; i<fMatchedClusterIndex->GetSize(); i++)
2319 if (fMatchedClusterIndex->At(i)==clsIndex) {
2320 Float_t r = TMath::Sqrt(fResidualEta->At(i)*fResidualEta->At(i) + fResidualPhi->At(i)*fResidualPhi->At(i));
2324 AliDebug(3,Form("Matched cluster index: index: %d, dEta: %2.4f, dPhi: %2.4f.\n",
2325 fMatchedClusterIndex->At(i),fResidualEta->At(i),fResidualPhi->At(i)));
2332 //____________________________________________________________________
2333 UInt_t AliEMCALRecoUtils::FindMatchedPosForTrack(Int_t trkIndex) const
2335 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2336 //Returns the position of the match in the fMatchedTrackIndex array
2337 Float_t tmpR = fCutR;
2340 for (Int_t i=0; i<fMatchedTrackIndex->GetSize(); i++)
2342 if (fMatchedTrackIndex->At(i)==trkIndex) {
2343 Float_t r = TMath::Sqrt(fResidualEta->At(i)*fResidualEta->At(i) + fResidualPhi->At(i)*fResidualPhi->At(i));
2347 AliDebug(3,Form("Matched track index: index: %d, dEta: %2.4f, dPhi: %2.4f.\n",
2348 fMatchedTrackIndex->At(i),fResidualEta->At(i),fResidualPhi->At(i)));
2355 //__________________________________________________________________________
2356 Bool_t AliEMCALRecoUtils::IsGoodCluster(AliVCluster *cluster,
2357 const AliEMCALGeometry *geom,
2358 AliVCaloCells* cells, Int_t bc)
2360 // check if the cluster survives some quality cut
2363 Bool_t isGood=kTRUE;
2365 if (!cluster || !cluster->IsEMCAL()) return kFALSE;
2366 if (ClusterContainsBadChannel(geom,cluster->GetCellsAbsId(),cluster->GetNCells())) return kFALSE;
2367 if (!CheckCellFiducialRegion(geom,cluster,cells)) return kFALSE;
2368 if (IsExoticCluster(cluster, cells,bc)) return kFALSE;
2373 //__________________________________________________________
2374 Bool_t AliEMCALRecoUtils::IsAccepted(AliESDtrack *esdTrack)
2376 // Given a esd track, return whether the track survive all the cuts
2378 // The different quality parameter are first
2379 // retrieved from the track. then it is found out what cuts the
2380 // track did not survive and finally the cuts are imposed.
2382 UInt_t status = esdTrack->GetStatus();
2384 Int_t nClustersITS = esdTrack->GetITSclusters(0);
2385 Int_t nClustersTPC = esdTrack->GetTPCclusters(0);
2387 Float_t chi2PerClusterITS = -1;
2388 Float_t chi2PerClusterTPC = -1;
2389 if (nClustersITS!=0)
2390 chi2PerClusterITS = esdTrack->GetITSchi2()/Float_t(nClustersITS);
2391 if (nClustersTPC!=0)
2392 chi2PerClusterTPC = esdTrack->GetTPCchi2()/Float_t(nClustersTPC);
2396 if (fTrackCutsType==kGlobalCut) {
2397 Float_t maxDCAToVertexXYPtDep = 0.0182 + 0.0350/TMath::Power(esdTrack->Pt(),1.01); //This expression comes from AliESDtrackCuts::GetStandardITSTPCTrackCuts2010()
2398 //AliDebug(3,Form("Track pT = %f, DCAtoVertexXY = %f",esdTrack->Pt(),MaxDCAToVertexXYPtDep));
2399 SetMaxDCAToVertexXY(maxDCAToVertexXYPtDep); //Set pT dependent DCA cut to vertex in x-y plane
2404 esdTrack->GetImpactParameters(b,bCov);
2405 if (bCov[0]<=0 || bCov[2]<=0) {
2406 AliDebug(1, "Estimated b resolution lower or equal zero!");
2407 bCov[0]=0; bCov[2]=0;
2410 Float_t dcaToVertexXY = b[0];
2411 Float_t dcaToVertexZ = b[1];
2412 Float_t dcaToVertex = -1;
2414 if (fCutDCAToVertex2D)
2415 dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY/fCutMaxDCAToVertexXY/fCutMaxDCAToVertexXY + dcaToVertexZ*dcaToVertexZ/fCutMaxDCAToVertexZ/fCutMaxDCAToVertexZ);
2417 dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY + dcaToVertexZ*dcaToVertexZ);
2421 Bool_t cuts[kNCuts];
2422 for (Int_t i=0; i<kNCuts; i++) cuts[i]=kFALSE;
2424 // track quality cuts
2425 if (fCutRequireTPCRefit && (status&AliESDtrack::kTPCrefit)==0)
2427 if (fCutRequireITSRefit && (status&AliESDtrack::kITSrefit)==0)
2429 if (nClustersTPC<fCutMinNClusterTPC)
2431 if (nClustersITS<fCutMinNClusterITS)
2433 if (chi2PerClusterTPC>fCutMaxChi2PerClusterTPC)
2435 if (chi2PerClusterITS>fCutMaxChi2PerClusterITS)
2437 if (!fCutAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0)
2439 if (fCutDCAToVertex2D && dcaToVertex > 1)
2441 if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexXY) > fCutMaxDCAToVertexXY)
2443 if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexZ) > fCutMaxDCAToVertexZ)
2446 if (fTrackCutsType==kGlobalCut) {
2447 //Require at least one SPD point + anything else in ITS
2448 if ( (esdTrack->HasPointOnITSLayer(0) || esdTrack->HasPointOnITSLayer(1)) == kFALSE)
2453 if (fCutRequireITSStandAlone || fCutRequireITSpureSA) {
2454 if ((status & AliESDtrack::kITSin) == 0 || (status & AliESDtrack::kTPCin)) {
2458 // ITS standalone tracks
2459 if (fCutRequireITSStandAlone && !fCutRequireITSpureSA) {
2460 if (status & AliESDtrack::kITSpureSA) cuts[11] = kTRUE;
2461 } else if (fCutRequireITSpureSA) {
2462 if (!(status & AliESDtrack::kITSpureSA)) cuts[11] = kTRUE;
2468 for (Int_t i=0; i<kNCuts; i++)
2469 if (cuts[i]) { cut = kTRUE ; }
2478 //_____________________________________
2479 void AliEMCALRecoUtils::InitTrackCuts()
2481 //Intilize the track cut criteria
2482 //By default these cuts are set according to AliESDtrackCuts::GetStandardTPCOnlyTrackCuts()
2483 //Also you can customize the cuts using the setters
2485 switch (fTrackCutsType)
2489 AliInfo(Form("Track cuts for matching: GetStandardTPCOnlyTrackCuts()"));
2491 SetMinNClustersTPC(70);
2492 SetMaxChi2PerClusterTPC(4);
2493 SetAcceptKinkDaughters(kFALSE);
2494 SetRequireTPCRefit(kFALSE);
2497 SetRequireITSRefit(kFALSE);
2498 SetMaxDCAToVertexZ(3.2);
2499 SetMaxDCAToVertexXY(2.4);
2500 SetDCAToVertex2D(kTRUE);
2507 AliInfo(Form("Track cuts for matching: GetStandardITSTPCTrackCuts2010(kTURE)"));
2509 SetMinNClustersTPC(70);
2510 SetMaxChi2PerClusterTPC(4);
2511 SetAcceptKinkDaughters(kFALSE);
2512 SetRequireTPCRefit(kTRUE);
2515 SetRequireITSRefit(kTRUE);
2516 SetMaxDCAToVertexZ(2);
2517 SetMaxDCAToVertexXY();
2518 SetDCAToVertex2D(kFALSE);
2525 AliInfo(Form("Track cuts for matching: Loose cut w/o DCA cut"));
2526 SetMinNClustersTPC(50);
2527 SetAcceptKinkDaughters(kTRUE);
2532 case kITSStandAlone:
2534 AliInfo(Form("Track cuts for matching: ITS Stand Alone tracks cut w/o DCA cut"));
2535 SetRequireITSRefit(kTRUE);
2536 SetRequireITSStandAlone(kTRUE);
2537 SetITSTrackSA(kTRUE);
2545 //________________________________________________________________________
2546 void AliEMCALRecoUtils::SetClusterMatchedToTrack(const AliVEvent *event)
2548 // Checks if tracks are matched to EMC clusters and set the matched EMCAL cluster index to ESD track.
2550 Int_t nTracks = event->GetNumberOfTracks();
2551 for (Int_t iTrack = 0; iTrack < nTracks; ++iTrack)
2553 AliVTrack* track = dynamic_cast<AliVTrack*>(event->GetTrack(iTrack));
2556 AliWarning(Form("Could not receive track %d", iTrack));
2560 Int_t matchClusIndex = GetMatchedClusterIndex(iTrack);
2561 track->SetEMCALcluster(matchClusIndex); //sets -1 if track not matched within residual
2562 /*the following can be done better if AliVTrack::SetStatus will be there. Patch pending with Andreas/Peter*/
2563 AliESDtrack* esdtrack = dynamic_cast<AliESDtrack*>(track);
2565 if (matchClusIndex != -1)
2566 esdtrack->SetStatus(AliESDtrack::kEMCALmatch);
2568 esdtrack->ResetStatus(AliESDtrack::kEMCALmatch);
2570 AliAODTrack* aodtrack = dynamic_cast<AliAODTrack*>(track);
2571 if (matchClusIndex != -1)
2572 aodtrack->SetStatus(AliESDtrack::kEMCALmatch);
2574 aodtrack->ResetStatus(AliESDtrack::kEMCALmatch);
2577 AliDebug(2,"Track matched to closest cluster");
2580 //_________________________________________________________________________
2581 void AliEMCALRecoUtils::SetTracksMatchedToCluster(const AliVEvent *event)
2583 // Checks if EMC clusters are matched to ESD track.
2584 // Adds track indexes of all the tracks matched to a cluster withing residuals in ESDCalocluster.
2586 for (Int_t iClus=0; iClus < event->GetNumberOfCaloClusters(); ++iClus)
2588 AliVCluster *cluster = event->GetCaloCluster(iClus);
2589 if (!cluster->IsEMCAL())
2592 Int_t nTracks = event->GetNumberOfTracks();
2593 TArrayI arrayTrackMatched(nTracks);
2595 // Get the closest track matched to the cluster
2597 Int_t matchTrackIndex = GetMatchedTrackIndex(iClus);
2598 if (matchTrackIndex != -1)
2600 arrayTrackMatched[nMatched] = matchTrackIndex;
2604 // Get all other tracks matched to the cluster
2605 for (Int_t iTrk=0; iTrk<nTracks; ++iTrk)
2607 AliVTrack* track = dynamic_cast<AliVTrack*>(event->GetTrack(iTrk));
2609 if( !track ) continue;
2611 if ( iTrk == matchTrackIndex ) continue;
2613 if ( track->GetEMCALcluster() == iClus )
2615 arrayTrackMatched[nMatched] = iTrk;
2620 //printf("Tender::SetTracksMatchedToCluster - cluster E %f, N matches %d, first match %d\n",cluster->E(),nMatched,arrayTrackMatched[0]);
2622 arrayTrackMatched.Set(nMatched);
2623 AliESDCaloCluster *esdcluster = dynamic_cast<AliESDCaloCluster*>(cluster);
2625 esdcluster->AddTracksMatched(arrayTrackMatched);
2626 else if (nMatched>0) {
2627 AliAODCaloCluster *aodcluster = dynamic_cast<AliAODCaloCluster*>(cluster);
2629 aodcluster->AddTrackMatched(event->GetTrack(arrayTrackMatched.At(0)));
2632 Float_t eta= -999, phi = -999;
2633 if (matchTrackIndex != -1)
2634 GetMatchedResiduals(iClus, eta, phi);
2635 cluster->SetTrackDistance(phi, eta);
2638 AliDebug(2,"Cluster matched to tracks");
2641 //___________________________________________________
2642 void AliEMCALRecoUtils::Print(const Option_t *) const
2646 printf("AliEMCALRecoUtils Settings: \n");
2647 printf("Misalignment shifts\n");
2648 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,
2649 fMisalTransShift[i*3],fMisalTransShift[i*3+1],fMisalTransShift[i*3+2],
2650 fMisalRotShift[i*3], fMisalRotShift[i*3+1], fMisalRotShift[i*3+2] );
2651 printf("Non linearity function %d, parameters:\n", fNonLinearityFunction);
2652 for (Int_t i=0; i<6; i++) printf("param[%d]=%f\n",i, fNonLinearityParams[i]);
2654 printf("Position Recalculation option %d, Particle Type %d, fW0 %2.2f, Recalibrate Data %d \n",fPosAlgo,fParticleType,fW0, fRecalibration);
2656 printf("Matching criteria: ");
2657 if (fCutEtaPhiSum) {
2658 printf("sqrt(dEta^2+dPhi^2)<%4.3f\n",fCutR);
2659 } else if (fCutEtaPhiSeparate) {
2660 printf("dEta<%4.3f, dPhi<%4.3f\n",fCutEta,fCutPhi);
2663 printf("please specify your cut criteria\n");
2664 printf("To cut on sqrt(dEta^2+dPhi^2), use: SwitchOnCutEtaPhiSum()\n");
2665 printf("To cut on dEta and dPhi separately, use: SwitchOnCutEtaPhiSeparate()\n");
2668 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);
2669 printf("Cluster selection window: dR < %2.0f\n",fClusterWindow);
2671 printf("Track cuts: \n");
2672 printf("Minimum track pT: %1.2f\n",fCutMinTrackPt);
2673 printf("AOD track selection: tpc only %d, or hybrid %d, or mask: %d\n",fAODTPCOnlyTracks,fAODHybridTracks, fAODFilterMask);
2674 printf("TPCRefit = %d, ITSRefit = %d\n",fCutRequireTPCRefit,fCutRequireITSRefit);
2675 printf("AcceptKinks = %d\n",fCutAcceptKinkDaughters);
2676 printf("MinNCulsterTPC = %d, MinNClusterITS = %d\n",fCutMinNClusterTPC,fCutMinNClusterITS);
2677 printf("MaxChi2TPC = %2.2f, MaxChi2ITS = %2.2f\n",fCutMaxChi2PerClusterTPC,fCutMaxChi2PerClusterITS);
2678 printf("DCSToVertex2D = %d, MaxDCAToVertexXY = %2.2f, MaxDCAToVertexZ = %2.2f\n",fCutDCAToVertex2D,fCutMaxDCAToVertexXY,fCutMaxDCAToVertexZ);