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(430.),
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(430.),
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;
226 if(reco.fResidualEta)
228 // assign or copy construct
231 *fResidualEta = *reco.fResidualEta;
235 fResidualEta = new TArrayF(*reco.fResidualEta);
240 if(fResidualEta)delete fResidualEta;
244 if(reco.fResidualPhi)
246 // assign or copy construct
249 *fResidualPhi = *reco.fResidualPhi;
253 fResidualPhi = new TArrayF(*reco.fResidualPhi);
258 if(fResidualPhi)delete fResidualPhi;
262 if(reco.fMatchedTrackIndex)
264 // assign or copy construct
265 if(fMatchedTrackIndex)
267 *fMatchedTrackIndex = *reco.fMatchedTrackIndex;
271 fMatchedTrackIndex = new TArrayI(*reco.fMatchedTrackIndex);
276 if(fMatchedTrackIndex)delete fMatchedTrackIndex;
277 fMatchedTrackIndex = 0;
280 if(reco.fMatchedClusterIndex)
282 // assign or copy construct
283 if(fMatchedClusterIndex)
285 *fMatchedClusterIndex = *reco.fMatchedClusterIndex;
289 fMatchedClusterIndex = new TArrayI(*reco.fMatchedClusterIndex);
294 if(fMatchedClusterIndex)delete fMatchedClusterIndex;
295 fMatchedClusterIndex = 0;
302 //_____________________________________
303 AliEMCALRecoUtils::~AliEMCALRecoUtils()
307 if(fEMCALRecalibrationFactors)
309 fEMCALRecalibrationFactors->Clear();
310 delete fEMCALRecalibrationFactors;
313 if(fEMCALTimeRecalibrationFactors)
315 fEMCALTimeRecalibrationFactors->Clear();
316 delete fEMCALTimeRecalibrationFactors;
319 if(fEMCALBadChannelMap)
321 fEMCALBadChannelMap->Clear();
322 delete fEMCALBadChannelMap;
325 delete fMatchedTrackIndex ;
326 delete fMatchedClusterIndex ;
327 delete fResidualEta ;
328 delete fResidualPhi ;
334 //_______________________________________________________________________________
335 Bool_t AliEMCALRecoUtils::AcceptCalibrateCell(const Int_t absID, const Int_t bc,
336 Float_t & amp, Double_t & time,
337 AliVCaloCells* cells)
339 // Reject cell if criteria not passed and calibrate it
341 AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
343 if(absID < 0 || absID >= 24*48*geom->GetNumberOfSuperModules()) return kFALSE;
345 Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1;
347 if(!geom->GetCellIndex(absID,imod,iTower,iIphi,iIeta))
349 // cell absID does not exist
354 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);
356 // Do not include bad channels found in analysis,
357 if( IsBadChannelsRemovalSwitchedOn() && GetEMCALChannelStatus(imod, ieta, iphi))
363 amp = cells->GetCellAmplitude(absID);
364 if(!fCellsRecalibrated && IsRecalibrationOn())
365 amp *= GetEMCALChannelRecalibrationFactor(imod,ieta,iphi);
369 time = cells->GetCellTime(absID);
371 RecalibrateCellTime(absID,bc,time);
376 //_____________________________________________________________________________
377 Bool_t AliEMCALRecoUtils::CheckCellFiducialRegion(const AliEMCALGeometry* geom,
378 const AliVCluster* cluster,
379 AliVCaloCells* cells)
381 // Given the list of AbsId of the cluster, get the maximum cell and
382 // check if there are fNCellsFromBorder from the calorimeter border
386 AliInfo("Cluster pointer null!");
390 //If the distance to the border is 0 or negative just exit accept all clusters
391 if(cells->GetType()==AliVCaloCells::kEMCALCell && fNCellsFromEMCALBorder <= 0 ) return kTRUE;
393 Int_t absIdMax = -1, iSM =-1, ieta = -1, iphi = -1;
394 Bool_t shared = kFALSE;
395 GetMaxEnergyCell(geom, cells, cluster, absIdMax, iSM, ieta, iphi, shared);
397 AliDebug(2,Form("Cluster Max AbsId %d, Cell Energy %2.2f, Cluster Energy %2.2f, Ncells from border %d, EMCAL eta=0 %d\n",
398 absIdMax, cells->GetCellAmplitude(absIdMax), cluster->E(), fNCellsFromEMCALBorder, fNoEMCALBorderAtEta0));
400 if(absIdMax==-1) return kFALSE;
402 //Check if the cell is close to the borders:
403 Bool_t okrow = kFALSE;
404 Bool_t okcol = kFALSE;
406 if(iSM < 0 || iphi < 0 || ieta < 0 )
408 AliFatal(Form("Negative value for super module: %d, or cell ieta: %d, or cell iphi: %d, check EMCAL geometry name\n",
415 if(iphi >= fNCellsFromEMCALBorder && iphi < 24-fNCellsFromEMCALBorder) okrow =kTRUE;
417 else if (iSM >=10 && ( ( geom->GetEMCGeometry()->GetGeoName()).Contains("12SMV1")))
419 if(iphi >= fNCellsFromEMCALBorder && iphi < 8-fNCellsFromEMCALBorder) okrow =kTRUE; //1/3 sm case
423 if(iphi >= fNCellsFromEMCALBorder && iphi < 12-fNCellsFromEMCALBorder) okrow =kTRUE; // half SM case
427 if(!fNoEMCALBorderAtEta0)
429 if(ieta > fNCellsFromEMCALBorder && ieta < 48-fNCellsFromEMCALBorder) okcol =kTRUE;
435 if(ieta >= fNCellsFromEMCALBorder) okcol = kTRUE;
439 if(ieta < 48-fNCellsFromEMCALBorder) okcol = kTRUE;
443 AliDebug(2,Form("EMCAL Cluster in %d cells fiducial volume: ieta %d, iphi %d, SM %d: column? %d, row? %d\nq",
444 fNCellsFromEMCALBorder, ieta, iphi, iSM, okcol, okrow));
448 //printf("Accept\n");
453 //printf("Reject\n");
454 AliDebug(2,Form("Reject cluster in border, max cell : ieta %d, iphi %d, SM %d\n",ieta, iphi, iSM));
461 //_______________________________________________________________________________
462 Bool_t AliEMCALRecoUtils::ClusterContainsBadChannel(const AliEMCALGeometry* geom,
463 const UShort_t* cellList,
466 // Check that in the cluster cells, there is no bad channel of those stored
467 // in fEMCALBadChannelMap or fPHOSBadChannelMap
469 if(!fRemoveBadChannels) return kFALSE;
470 if(!fEMCALBadChannelMap) return kFALSE;
475 for(Int_t iCell = 0; iCell<nCells; iCell++)
477 //Get the column and row
478 Int_t iTower = -1, iIphi = -1, iIeta = -1;
479 geom->GetCellIndex(cellList[iCell],imod,iTower,iIphi,iIeta);
480 if(fEMCALBadChannelMap->GetEntries() <= imod) continue;
481 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,irow,icol);
482 if(GetEMCALChannelStatus(imod, icol, irow))
484 AliDebug(2,Form("Cluster with bad channel: SM %d, col %d, row %d\n",imod, icol, irow));
488 }// cell cluster loop
494 //___________________________________________________________________________
495 Float_t AliEMCALRecoUtils::GetECross(const Int_t absID, const Double_t tcell,
496 AliVCaloCells* cells, const Int_t bc)
498 //Calculate the energy in the cross around the energy given cell
500 AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
502 Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1;
503 geom->GetCellIndex(absID,imod,iTower,iIphi,iIeta);
504 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);
506 //Get close cells index, energy and time, not in corners
511 if( iphi < AliEMCALGeoParams::fgkEMCALRows-1) absID1 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi+1, ieta);
512 if( iphi > 0 ) absID2 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi-1, ieta);
514 // In case of cell in eta = 0 border, depending on SM shift the cross cell index
519 if ( ieta == AliEMCALGeoParams::fgkEMCALCols-1 && !(imod%2) )
521 absID3 = geom-> GetAbsCellIdFromCellIndexes(imod+1, iphi, 0);
522 absID4 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta-1);
524 else if( ieta == 0 && imod%2 )
526 absID3 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta+1);
527 absID4 = geom-> GetAbsCellIdFromCellIndexes(imod-1, iphi, AliEMCALGeoParams::fgkEMCALCols-1);
531 if( ieta < AliEMCALGeoParams::fgkEMCALCols-1 )
532 absID3 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta+1);
534 absID4 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta-1);
537 //printf("IMOD %d, AbsId %d, a %d, b %d, c %d e %d \n",imod,absID,absID1,absID2,absID3,absID4);
539 Float_t ecell1 = 0, ecell2 = 0, ecell3 = 0, ecell4 = 0;
540 Double_t tcell1 = 0, tcell2 = 0, tcell3 = 0, tcell4 = 0;
541 Bool_t accept1 = 0, accept2 = 0, accept3 = 0, accept4 = 0;
543 accept1 = AcceptCalibrateCell(absID1,bc, ecell1,tcell1,cells);
544 accept2 = AcceptCalibrateCell(absID2,bc, ecell2,tcell2,cells);
545 accept3 = AcceptCalibrateCell(absID3,bc, ecell3,tcell3,cells);
546 accept4 = AcceptCalibrateCell(absID4,bc, ecell4,tcell4,cells);
549 printf("Cell absID %d \n",absID);
550 printf("\t accept1 %d, accept2 %d, accept3 %d, accept4 %d\n",
551 accept1,accept2,accept3,accept4);
552 printf("\t id %d: id1 %d, id2 %d, id3 %d, id4 %d\n",
553 absID,absID1,absID2,absID3,absID4);
554 printf("\t e %f: e1 %f, e2 %f, e3 %f, e4 %f\n",
555 ecell,ecell1,ecell2,ecell3,ecell4);
556 printf("\t t %f: t1 %f, t2 %f, t3 %f, t4 %f;\n dt1 %f, dt2 %f, dt3 %f, dt4 %f\n",
557 tcell*1.e9,tcell1*1.e9,tcell2*1.e9,tcell3*1.e9,tcell4*1.e9,
558 TMath::Abs(tcell-tcell1)*1.e9, TMath::Abs(tcell-tcell2)*1.e9, TMath::Abs(tcell-tcell3)*1.e9, TMath::Abs(tcell-tcell4)*1.e9);
561 if(TMath::Abs(tcell-tcell1)*1.e9 > fExoticCellDiffTime) ecell1 = 0 ;
562 if(TMath::Abs(tcell-tcell2)*1.e9 > fExoticCellDiffTime) ecell2 = 0 ;
563 if(TMath::Abs(tcell-tcell3)*1.e9 > fExoticCellDiffTime) ecell3 = 0 ;
564 if(TMath::Abs(tcell-tcell4)*1.e9 > fExoticCellDiffTime) ecell4 = 0 ;
566 return ecell1+ecell2+ecell3+ecell4;
570 //_____________________________________________________________________________________________
571 Bool_t AliEMCALRecoUtils::IsExoticCell(const Int_t absID, AliVCaloCells* cells, const Int_t bc)
573 // Look to cell neighbourhood and reject if it seems exotic
574 // Do before recalibrating the cells
576 if(!fRejectExoticCells) return kFALSE;
580 Bool_t accept = AcceptCalibrateCell(absID, bc, ecell ,tcell ,cells);
582 if(!accept) return kTRUE; // reject this cell
584 if(ecell < fExoticCellMinAmplitude) return kFALSE; // do not reject low energy cells
586 Float_t eCross = GetECross(absID,tcell,cells,bc);
588 if(1-eCross/ecell > fExoticCellFraction)
590 AliDebug(2,Form("AliEMCALRecoUtils::IsExoticCell() - EXOTIC CELL id %d, eCell %f, eCross %f, 1-eCross/eCell %f\n",
591 absID,ecell,eCross,1-eCross/ecell));
598 //___________________________________________________________________
599 Bool_t AliEMCALRecoUtils::IsExoticCluster(const AliVCluster *cluster,
600 AliVCaloCells *cells,
603 // Check if the cluster highest energy tower is exotic
607 AliInfo("Cluster pointer null!");
611 if(!fRejectExoticCluster) return kFALSE;
613 // Get highest energy tower
614 AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
615 Int_t iSupMod = -1, absId = -1, ieta = -1, iphi = -1;
616 Bool_t shared = kFALSE;
617 GetMaxEnergyCell(geom, cells, cluster, absId, iSupMod, ieta, iphi, shared);
619 return IsExoticCell(absId,cells,bc);
623 //_______________________________________________________________________
624 Float_t AliEMCALRecoUtils::SmearClusterEnergy(const AliVCluster* cluster)
626 //In case of MC analysis, smear energy to match resolution/calibration in real data
630 AliInfo("Cluster pointer null!");
634 Float_t energy = cluster->E() ;
635 Float_t rdmEnergy = energy ;
636 if(fSmearClusterEnergy)
638 rdmEnergy = fRandom.Gaus(energy,fSmearClusterParam[0] * TMath::Sqrt(energy) +
639 fSmearClusterParam[1] * energy +
640 fSmearClusterParam[2] );
641 AliDebug(2, Form("Energy: original %f, smeared %f\n", energy, rdmEnergy));
647 //____________________________________________________________________________
648 Float_t AliEMCALRecoUtils::CorrectClusterEnergyLinearity(AliVCluster* cluster)
650 // Correct cluster energy from non linearity functions
654 AliInfo("Cluster pointer null!");
658 Float_t energy = cluster->E();
662 // Clusters with less than 50 MeV or negative are not possible
663 AliInfo(Form("Too Low Cluster energy!, E = %f < 0.05 GeV",energy));
667 switch (fNonLinearityFunction)
672 //Non-Linearity correction (from MC with function ([0]*exp(-[1]/E))+(([2]/([3]*2.*TMath::Pi())*exp(-(E-[4])^2/(2.*[3]^2)))))
673 //fNonLinearityParams[0] = 1.014;
674 //fNonLinearityParams[1] =-0.03329;
675 //fNonLinearityParams[2] =-0.3853;
676 //fNonLinearityParams[3] = 0.5423;
677 //fNonLinearityParams[4] =-0.4335;
678 energy *= (fNonLinearityParams[0]*exp(-fNonLinearityParams[1]/energy))+
679 ((fNonLinearityParams[2]/(fNonLinearityParams[3]*2.*TMath::Pi())*
680 exp(-(energy-fNonLinearityParams[4])*(energy-fNonLinearityParams[4])/(2.*fNonLinearityParams[3]*fNonLinearityParams[3]))));
686 //Non-Linearity correction (from MC with function [0]/((x+[1])^[2]))+1;
687 //fNonLinearityParams[0] = 3.11111e-02;
688 //fNonLinearityParams[1] =-5.71666e-02;
689 //fNonLinearityParams[2] = 5.67995e-01;
691 energy *= fNonLinearityParams[0]/TMath::Power(energy+fNonLinearityParams[1],fNonLinearityParams[2])+1;
697 //Same as beam test corrected, change parameters
698 //fNonLinearityParams[0] = 9.81039e-01
699 //fNonLinearityParams[1] = 1.13508e-01;
700 //fNonLinearityParams[2] = 1.00173e+00;
701 //fNonLinearityParams[3] = 9.67998e-02;
702 //fNonLinearityParams[4] = 2.19381e+02;
703 //fNonLinearityParams[5] = 6.31604e+01;
704 //fNonLinearityParams[6] = 1;
705 energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
713 //Non-Linearity correction (from Olga Data with function p0+p1*exp(-p2*E))
714 //fNonLinearityParams[0] = 1.04;
715 //fNonLinearityParams[1] = -0.1445;
716 //fNonLinearityParams[2] = 1.046;
717 energy /= (fNonLinearityParams[0]+fNonLinearityParams[1]*exp(-fNonLinearityParams[2]*energy)); //Olga function
721 case kPi0GammaConversion:
723 //Non-Linearity correction (Nicolas from Dimitri Data with function C*[1-a*exp(-b*E)])
724 //fNonLinearityParams[0] = 0.139393/0.1349766;
725 //fNonLinearityParams[1] = 0.0566186;
726 //fNonLinearityParams[2] = 0.982133;
727 energy /= fNonLinearityParams[0]*(1-fNonLinearityParams[1]*exp(-fNonLinearityParams[2]*energy));
734 //From beam test, Alexei's results, for different ZS thresholds
735 // th=30 MeV; th = 45 MeV; th = 75 MeV
736 //fNonLinearityParams[0] = 1.007; 1.003; 1.002
737 //fNonLinearityParams[1] = 0.894; 0.719; 0.797
738 //fNonLinearityParams[2] = 0.246; 0.334; 0.358
739 //Rescale the param[0] with 1.03
740 energy /= fNonLinearityParams[0]/(1+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]));
745 case kBeamTestCorrected:
747 //From beam test, corrected for material between beam and EMCAL
748 //fNonLinearityParams[0] = 0.99078
749 //fNonLinearityParams[1] = 0.161499;
750 //fNonLinearityParams[2] = 0.655166;
751 //fNonLinearityParams[3] = 0.134101;
752 //fNonLinearityParams[4] = 163.282;
753 //fNonLinearityParams[5] = 23.6904;
754 //fNonLinearityParams[6] = 0.978;
755 energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
760 case kBeamTestCorrectedv2:
762 //From beam test, corrected for material between beam and EMCAL
763 //fNonLinearityParams[0] = 0.983504;
764 //fNonLinearityParams[1] = 0.210106;
765 //fNonLinearityParams[2] = 0.897274;
766 //fNonLinearityParams[3] = 0.0829064;
767 //fNonLinearityParams[4] = 152.299;
768 //fNonLinearityParams[5] = 31.5028;
769 //fNonLinearityParams[6] = 0.968;
770 energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
776 AliDebug(2,"No correction on the energy\n");
784 //__________________________________________________
785 void AliEMCALRecoUtils::InitNonLinearityParam()
787 //Initialising Non Linearity Parameters
789 if(fNonLinearityFunction == kPi0MC)
791 fNonLinearityParams[0] = 1.014;
792 fNonLinearityParams[1] = -0.03329;
793 fNonLinearityParams[2] = -0.3853;
794 fNonLinearityParams[3] = 0.5423;
795 fNonLinearityParams[4] = -0.4335;
798 if(fNonLinearityFunction == kPi0MCv2)
800 fNonLinearityParams[0] = 3.11111e-02;
801 fNonLinearityParams[1] =-5.71666e-02;
802 fNonLinearityParams[2] = 5.67995e-01;
805 if(fNonLinearityFunction == kPi0MCv3)
807 fNonLinearityParams[0] = 9.81039e-01;
808 fNonLinearityParams[1] = 1.13508e-01;
809 fNonLinearityParams[2] = 1.00173e+00;
810 fNonLinearityParams[3] = 9.67998e-02;
811 fNonLinearityParams[4] = 2.19381e+02;
812 fNonLinearityParams[5] = 6.31604e+01;
813 fNonLinearityParams[6] = 1;
816 if(fNonLinearityFunction == kPi0GammaGamma)
818 fNonLinearityParams[0] = 1.04;
819 fNonLinearityParams[1] = -0.1445;
820 fNonLinearityParams[2] = 1.046;
823 if(fNonLinearityFunction == kPi0GammaConversion)
825 fNonLinearityParams[0] = 0.139393;
826 fNonLinearityParams[1] = 0.0566186;
827 fNonLinearityParams[2] = 0.982133;
830 if(fNonLinearityFunction == kBeamTest)
832 if(fNonLinearThreshold == 30)
834 fNonLinearityParams[0] = 1.007;
835 fNonLinearityParams[1] = 0.894;
836 fNonLinearityParams[2] = 0.246;
838 if(fNonLinearThreshold == 45)
840 fNonLinearityParams[0] = 1.003;
841 fNonLinearityParams[1] = 0.719;
842 fNonLinearityParams[2] = 0.334;
844 if(fNonLinearThreshold == 75)
846 fNonLinearityParams[0] = 1.002;
847 fNonLinearityParams[1] = 0.797;
848 fNonLinearityParams[2] = 0.358;
852 if(fNonLinearityFunction == kBeamTestCorrected)
854 fNonLinearityParams[0] = 0.99078;
855 fNonLinearityParams[1] = 0.161499;
856 fNonLinearityParams[2] = 0.655166;
857 fNonLinearityParams[3] = 0.134101;
858 fNonLinearityParams[4] = 163.282;
859 fNonLinearityParams[5] = 23.6904;
860 fNonLinearityParams[6] = 0.978;
863 if(fNonLinearityFunction == kBeamTestCorrectedv2)
865 fNonLinearityParams[0] = 0.983504;
866 fNonLinearityParams[1] = 0.210106;
867 fNonLinearityParams[2] = 0.897274;
868 fNonLinearityParams[3] = 0.0829064;
869 fNonLinearityParams[4] = 152.299;
870 fNonLinearityParams[5] = 31.5028;
871 fNonLinearityParams[6] = 0.968;
875 //_________________________________________________________
876 Float_t AliEMCALRecoUtils::GetDepth(const Float_t energy,
877 const Int_t iParticle,
878 const Int_t iSM) const
880 //Calculate shower depth for a given cluster energy and particle type
886 Float_t arg = energy*1000/ ecr; //Multiply energy by 1000 to transform to MeV
894 depth = x0 * (TMath::Log(arg) + 0.5);
901 depth = x0 * (TMath::Log(arg) - 0.5);
909 gGeoManager->cd("ALIC_1/XEN1_1");
910 TGeoNode *geoXEn1 = gGeoManager->GetCurrentNode();
911 TGeoNodeMatrix *geoSM = dynamic_cast<TGeoNodeMatrix *>(geoXEn1->GetDaughter(iSM));
914 TGeoVolume *geoSMVol = geoSM->GetVolume();
915 TGeoShape *geoSMShape = geoSMVol->GetShape();
916 TGeoBBox *geoBox = dynamic_cast<TGeoBBox *>(geoSMShape);
917 if(geoBox) depth = 0.5 * geoBox->GetDX()*2 ;
918 else AliFatal("Null GEANT box");
920 else AliFatal("NULL GEANT node matrix");
927 depth = x0 * (TMath::Log(arg) - 0.5);
936 depth = x0 * (TMath::Log(arg) + 0.5);
942 //____________________________________________________________________
943 void AliEMCALRecoUtils::GetMaxEnergyCell(const AliEMCALGeometry *geom,
944 AliVCaloCells* cells,
945 const AliVCluster* clu,
952 //For a given CaloCluster gets the absId of the cell
953 //with maximum energy deposit.
956 Double_t eCell = -1.;
957 Float_t fraction = 1.;
958 Float_t recalFactor = 1.;
959 Int_t cellAbsId = -1 ;
968 AliInfo("Cluster pointer null!");
969 absId=-1; iSupMod0=-1, ieta = -1; iphi = -1; shared = -1;
973 for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++)
975 cellAbsId = clu->GetCellAbsId(iDig);
976 fraction = clu->GetCellAmplitudeFraction(iDig);
977 //printf("a Cell %d, id, %d, amp %f, fraction %f\n",iDig,cellAbsId,cells->GetCellAmplitude(cellAbsId),fraction);
978 if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
979 geom->GetCellIndex(cellAbsId,iSupMod,iTower,iIphi,iIeta);
980 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);
985 else if(iSupMod0!=iSupMod)
988 //printf("AliEMCALRecoUtils::GetMaxEnergyCell() - SHARED CLUSTER\n");
990 if(!fCellsRecalibrated && IsRecalibrationOn())
992 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
994 eCell = cells->GetCellAmplitude(cellAbsId)*fraction*recalFactor;
995 //printf("b Cell %d, id, %d, amp %f, fraction %f\n",iDig,cellAbsId,eCell,fraction);
1000 //printf("\t new max: cell %d, e %f, ecell %f\n",maxId, eMax,eCell);
1004 //Get from the absid the supermodule, tower and eta/phi numbers
1005 geom->GetCellIndex(absId,iSupMod,iTower,iIphi,iIeta);
1006 //Gives SuperModule and Tower numbers
1007 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,
1008 iIphi, iIeta,iphi,ieta);
1009 //printf("Max id %d, iSM %d, col %d, row %d\n",absId,iSupMod,ieta,iphi);
1010 //printf("Max end---\n");
1013 //______________________________________
1014 void AliEMCALRecoUtils::InitParameters()
1016 // Initialize data members with default values
1018 fParticleType = kPhoton;
1019 fPosAlgo = kUnchanged;
1022 fNonLinearityFunction = kNoCorrection;
1023 fNonLinearThreshold = 30;
1025 fExoticCellFraction = 0.97;
1026 fExoticCellDiffTime = 1e6;
1027 fExoticCellMinAmplitude = 0.5;
1029 fAODFilterMask = 128;
1030 fAODHybridTracks = kFALSE;
1031 fAODTPCOnlyTracks = kTRUE;
1033 fCutEtaPhiSum = kTRUE;
1034 fCutEtaPhiSeparate = kFALSE;
1040 fClusterWindow = 100;
1045 fTrackCutsType = kLooseCut;
1048 fCutMinNClusterTPC = -1;
1049 fCutMinNClusterITS = -1;
1051 fCutMaxChi2PerClusterTPC = 1e10;
1052 fCutMaxChi2PerClusterITS = 1e10;
1054 fCutRequireTPCRefit = kFALSE;
1055 fCutRequireITSRefit = kFALSE;
1056 fCutAcceptKinkDaughters = kFALSE;
1058 fCutMaxDCAToVertexXY = 1e10;
1059 fCutMaxDCAToVertexZ = 1e10;
1060 fCutDCAToVertex2D = kFALSE;
1062 fCutRequireITSStandAlone = kFALSE; //MARCEL
1063 fCutRequireITSpureSA = kFALSE; //Marcel
1065 //Misalignment matrices
1066 for(Int_t i = 0; i < 15 ; i++)
1068 fMisalTransShift[i] = 0.;
1069 fMisalRotShift[i] = 0.;
1073 for(Int_t i = 0; i < 7 ; i++) fNonLinearityParams[i] = 0.;
1075 //For kBeamTestCorrectedv2 case, but default is no correction
1076 fNonLinearityParams[0] = 0.983504;
1077 fNonLinearityParams[1] = 0.210106;
1078 fNonLinearityParams[2] = 0.897274;
1079 fNonLinearityParams[3] = 0.0829064;
1080 fNonLinearityParams[4] = 152.299;
1081 fNonLinearityParams[5] = 31.5028;
1082 fNonLinearityParams[6] = 0.968;
1084 //Cluster energy smearing
1085 fSmearClusterEnergy = kFALSE;
1086 fSmearClusterParam[0] = 0.07; // * sqrt E term
1087 fSmearClusterParam[1] = 0.00; // * E term
1088 fSmearClusterParam[2] = 0.00; // constant
1091 //_____________________________________________________
1092 void AliEMCALRecoUtils::InitEMCALRecalibrationFactors()
1094 //Init EMCAL recalibration factors
1095 AliDebug(2,"AliCalorimeterUtils::InitEMCALRecalibrationFactors()");
1096 //In order to avoid rewriting the same histograms
1097 Bool_t oldStatus = TH1::AddDirectoryStatus();
1098 TH1::AddDirectory(kFALSE);
1100 fEMCALRecalibrationFactors = new TObjArray(12);
1101 for (int i = 0; i < 12; i++)
1102 fEMCALRecalibrationFactors->Add(new TH2F(Form("EMCALRecalFactors_SM%d",i),
1103 Form("EMCALRecalFactors_SM%d",i), 48, 0, 48, 24, 0, 24));
1104 //Init the histograms with 1
1105 for (Int_t sm = 0; sm < 12; sm++)
1107 for (Int_t i = 0; i < 48; i++)
1109 for (Int_t j = 0; j < 24; j++)
1111 SetEMCALChannelRecalibrationFactor(sm,i,j,1.);
1116 fEMCALRecalibrationFactors->SetOwner(kTRUE);
1117 fEMCALRecalibrationFactors->Compress();
1119 //In order to avoid rewriting the same histograms
1120 TH1::AddDirectory(oldStatus);
1123 //_________________________________________________________
1124 void AliEMCALRecoUtils::InitEMCALTimeRecalibrationFactors()
1126 //Init EMCAL recalibration factors
1127 AliDebug(2,"AliCalorimeterUtils::InitEMCALRecalibrationFactors()");
1128 //In order to avoid rewriting the same histograms
1129 Bool_t oldStatus = TH1::AddDirectoryStatus();
1130 TH1::AddDirectory(kFALSE);
1132 fEMCALTimeRecalibrationFactors = new TObjArray(4);
1133 for (int i = 0; i < 4; i++)
1134 fEMCALTimeRecalibrationFactors->Add(new TH1F(Form("hAllTimeAvBC%d",i),
1135 Form("hAllTimeAvBC%d",i),
1136 48*24*12,0.,48*24*12) );
1137 //Init the histograms with 1
1138 for (Int_t bc = 0; bc < 4; bc++)
1140 for (Int_t i = 0; i < 48*24*12; i++)
1141 SetEMCALChannelTimeRecalibrationFactor(bc,i,0.);
1144 fEMCALTimeRecalibrationFactors->SetOwner(kTRUE);
1145 fEMCALTimeRecalibrationFactors->Compress();
1147 //In order to avoid rewriting the same histograms
1148 TH1::AddDirectory(oldStatus);
1151 //____________________________________________________
1152 void AliEMCALRecoUtils::InitEMCALBadChannelStatusMap()
1154 //Init EMCAL bad channels map
1155 AliDebug(2,"AliEMCALRecoUtils::InitEMCALBadChannelStatusMap()");
1156 //In order to avoid rewriting the same histograms
1157 Bool_t oldStatus = TH1::AddDirectoryStatus();
1158 TH1::AddDirectory(kFALSE);
1160 fEMCALBadChannelMap = new TObjArray(12);
1161 //TH2F * hTemp = new TH2I("EMCALBadChannelMap","EMCAL SuperModule bad channel map", 48, 0, 48, 24, 0, 24);
1162 for (int i = 0; i < 12; i++)
1164 fEMCALBadChannelMap->Add(new TH2I(Form("EMCALBadChannelMap_Mod%d",i),Form("EMCALBadChannelMap_Mod%d",i), 48, 0, 48, 24, 0, 24));
1167 fEMCALBadChannelMap->SetOwner(kTRUE);
1168 fEMCALBadChannelMap->Compress();
1170 //In order to avoid rewriting the same histograms
1171 TH1::AddDirectory(oldStatus);
1174 //____________________________________________________________________________
1175 void AliEMCALRecoUtils::RecalibrateClusterEnergy(const AliEMCALGeometry* geom,
1176 AliVCluster * cluster,
1177 AliVCaloCells * cells,
1180 // Recalibrate the cluster energy and Time, considering the recalibration map
1181 // and the energy of the cells and time that compose the cluster.
1182 // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1186 AliInfo("Cluster pointer null!");
1190 //Get the cluster number of cells and list of absId, check what kind of cluster do we have.
1191 UShort_t * index = cluster->GetCellsAbsId() ;
1192 Double_t * fraction = cluster->GetCellsAmplitudeFraction() ;
1193 Int_t ncells = cluster->GetNCells();
1195 //Initialize some used variables
1198 Int_t icol =-1, irow =-1, imod=1;
1199 Float_t factor = 1, frac = 0;
1200 Int_t absIdMax = -1;
1203 //Loop on the cells, get the cell amplitude and recalibration factor, multiply and and to the new energy
1204 for(Int_t icell = 0; icell < ncells; icell++)
1206 absId = index[icell];
1207 frac = fraction[icell];
1208 if(frac < 1e-5) frac = 1; //in case of EMCAL, this is set as 0 since unfolding is off
1210 if(!fCellsRecalibrated && IsRecalibrationOn())
1213 Int_t iTower = -1, iIphi = -1, iIeta = -1;
1214 geom->GetCellIndex(absId,imod,iTower,iIphi,iIeta);
1215 if(fEMCALRecalibrationFactors->GetEntries() <= imod) continue;
1216 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,irow,icol);
1217 factor = GetEMCALChannelRecalibrationFactor(imod,icol,irow);
1219 AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - recalibrate cell: module %d, col %d, row %d, cell fraction %f,recalibration factor %f, cell energy %f\n",
1220 imod,icol,irow,frac,factor,cells->GetCellAmplitude(absId)));
1224 energy += cells->GetCellAmplitude(absId)*factor*frac;
1226 if(emax < cells->GetCellAmplitude(absId)*factor*frac)
1228 emax = cells->GetCellAmplitude(absId)*factor*frac;
1233 AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - Energy before %f, after %f \n",cluster->E(),energy));
1235 cluster->SetE(energy);
1237 // Recalculate time of cluster
1238 Double_t timeorg = cluster->GetTOF();
1240 Double_t time = cells->GetCellTime(absIdMax);
1241 if(!fCellsRecalibrated && IsTimeRecalibrationOn())
1242 RecalibrateCellTime(absIdMax,bc,time);
1244 cluster->SetTOF(time);
1246 AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - Time before %f, after %f \n",timeorg,cluster->GetTOF()));
1249 //_____________________________________________________________
1250 void AliEMCALRecoUtils::RecalibrateCells(AliVCaloCells * cells,
1253 // Recalibrate the cells time and energy, considering the recalibration map and the energy
1254 // of the cells that compose the cluster.
1255 // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1257 if(!IsRecalibrationOn() && !IsTimeRecalibrationOn() && !IsBadChannelsRemovalSwitchedOn()) return;
1261 AliInfo("Cells pointer null!");
1266 Bool_t accept = kFALSE;
1269 Double_t ecellin = 0;
1270 Double_t tcellin = 0;
1274 Int_t nEMcell = cells->GetNumberOfCells() ;
1275 for (Int_t iCell = 0; iCell < nEMcell; iCell++)
1277 cells->GetCell( iCell, absId, ecellin, tcellin, mclabel, efrac );
1279 accept = AcceptCalibrateCell(absId, bc, ecell ,tcell ,cells);
1287 cells->SetCell(iCell,absId,ecell, tcell, mclabel, efrac);
1290 fCellsRecalibrated = kTRUE;
1293 //_______________________________________________________________________________________________________
1294 void AliEMCALRecoUtils::RecalibrateCellTime(const Int_t absId, const Int_t bc, Double_t & celltime) const
1296 // Recalibrate time of cell with absID considering the recalibration map
1297 // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1299 if(!fCellsRecalibrated && IsTimeRecalibrationOn() && bc >= 0)
1301 celltime -= GetEMCALChannelTimeRecalibrationFactor(bc%4,absId)*1.e-9; ;
1305 //______________________________________________________________________________
1306 void AliEMCALRecoUtils::RecalculateClusterPosition(const AliEMCALGeometry *geom,
1307 AliVCaloCells* cells,
1310 //For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1314 AliInfo("Cluster pointer null!");
1318 if (fPosAlgo==kPosTowerGlobal) RecalculateClusterPositionFromTowerGlobal( geom, cells, clu);
1319 else if(fPosAlgo==kPosTowerIndex) RecalculateClusterPositionFromTowerIndex ( geom, cells, clu);
1320 else AliDebug(2,"Algorithm to recalculate position not selected, do nothing.");
1323 //_____________________________________________________________________________________________
1324 void AliEMCALRecoUtils::RecalculateClusterPositionFromTowerGlobal(const AliEMCALGeometry *geom,
1325 AliVCaloCells* cells,
1328 // For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1329 // The algorithm is a copy of what is done in AliEMCALRecPoint
1331 Double_t eCell = 0.;
1332 Float_t fraction = 1.;
1333 Float_t recalFactor = 1.;
1336 Int_t iTower = -1, iIphi = -1, iIeta = -1;
1337 Int_t iSupModMax = -1, iSM=-1, iphi = -1, ieta = -1;
1338 Float_t weight = 0., totalWeight=0.;
1339 Float_t newPos[3] = {0,0,0};
1340 Double_t pLocal[3], pGlobal[3];
1341 Bool_t shared = kFALSE;
1343 Float_t clEnergy = clu->E(); //Energy already recalibrated previously
1346 GetMaxEnergyCell(geom, cells, clu, absId, iSupModMax, ieta, iphi,shared);
1347 Double_t depth = GetDepth(clEnergy,fParticleType,iSupModMax) ;
1349 //printf("** Cluster energy %f, ncells %d, depth %f\n",clEnergy,clu->GetNCells(),depth);
1351 for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++)
1353 absId = clu->GetCellAbsId(iDig);
1354 fraction = clu->GetCellAmplitudeFraction(iDig);
1355 if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1357 if (!fCellsRecalibrated)
1359 geom->GetCellIndex(absId,iSM,iTower,iIphi,iIeta);
1360 geom->GetCellPhiEtaIndexInSModule(iSM,iTower,iIphi, iIeta,iphi,ieta);
1362 if(IsRecalibrationOn())
1364 recalFactor = GetEMCALChannelRecalibrationFactor(iSM,ieta,iphi);
1368 eCell = cells->GetCellAmplitude(absId)*fraction*recalFactor;
1370 weight = GetCellWeight(eCell,clEnergy);
1371 totalWeight += weight;
1373 geom->RelPosCellInSModule(absId,depth,pLocal[0],pLocal[1],pLocal[2]);
1374 //printf("pLocal (%f,%f,%f), SM %d, absId %d\n",pLocal[0],pLocal[1],pLocal[2],iSupModMax,absId);
1375 geom->GetGlobal(pLocal,pGlobal,iSupModMax);
1376 //printf("pLocal (%f,%f,%f)\n",pGlobal[0],pGlobal[1],pGlobal[2]);
1378 for(int i=0; i<3; i++ ) newPos[i] += (weight*pGlobal[i]);
1383 for(int i=0; i<3; i++ ) newPos[i] /= totalWeight;
1386 //Float_t pos[]={0,0,0};
1387 //clu->GetPosition(pos);
1388 //printf("OldPos : %2.3f,%2.3f,%2.3f\n",pos[0],pos[1],pos[2]);
1389 //printf("NewPos : %2.3f,%2.3f,%2.3f\n",newPos[0],newPos[1],newPos[2]);
1391 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]);
1399 newPos[0] +=fMisalTransShift[0];//+=1.134;
1400 newPos[1] +=fMisalTransShift[1];//+=8.2;
1401 newPos[2] +=fMisalTransShift[2];//+=1.197;
1402 //printf(" + : %2.3f,%2.3f,%2.3f\n",fMisalTransShift[0],fMisalTransShift[1],fMisalTransShift[2]);
1404 //printf("NewPos : %2.3f,%2.3f,%2.3f\n",newPos[0],newPos[1],newPos[2]);
1406 clu->SetPosition(newPos);
1409 //____________________________________________________________________________________________
1410 void AliEMCALRecoUtils::RecalculateClusterPositionFromTowerIndex(const AliEMCALGeometry *geom,
1411 AliVCaloCells* cells,
1414 // For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1415 // The algorithm works with the tower indeces, averages the indeces and from them it calculates the global position
1417 Double_t eCell = 1.;
1418 Float_t fraction = 1.;
1419 Float_t recalFactor = 1.;
1423 Int_t iIphi = -1, iIeta = -1;
1424 Int_t iSupMod = -1, iSupModMax = -1;
1425 Int_t iphi = -1, ieta =-1;
1426 Bool_t shared = kFALSE;
1428 Float_t clEnergy = clu->E(); //Energy already recalibrated previously.
1432 GetMaxEnergyCell(geom, cells, clu, absId, iSupModMax, ieta, iphi,shared);
1433 Float_t depth = GetDepth(clEnergy,fParticleType,iSupMod) ;
1435 Float_t weight = 0., weightedCol = 0., weightedRow = 0., totalWeight=0.;
1436 Bool_t areInSameSM = kTRUE; //exclude clusters with cells in different SMs for now
1437 Int_t startingSM = -1;
1439 for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++)
1441 absId = clu->GetCellAbsId(iDig);
1442 fraction = clu->GetCellAmplitudeFraction(iDig);
1443 if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1445 if (iDig==0) startingSM = iSupMod;
1446 else if(iSupMod != startingSM) areInSameSM = kFALSE;
1448 eCell = cells->GetCellAmplitude(absId);
1450 geom->GetCellIndex(absId,iSupMod,iTower,iIphi,iIeta);
1451 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);
1453 if (!fCellsRecalibrated)
1455 if(IsRecalibrationOn())
1457 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1461 eCell = cells->GetCellAmplitude(absId)*fraction*recalFactor;
1463 weight = GetCellWeight(eCell,clEnergy);
1464 if(weight < 0) weight = 0;
1465 totalWeight += weight;
1466 weightedCol += ieta*weight;
1467 weightedRow += iphi*weight;
1469 //printf("Max cell? cell %d, amplitude org %f, fraction %f, recalibration %f, amplitude new %f \n",cellAbsId, cells->GetCellAmplitude(cellAbsId), fraction, recalFactor, eCell) ;
1472 Float_t xyzNew[]={0.,0.,0.};
1473 if(areInSameSM == kTRUE)
1475 //printf("In Same SM\n");
1476 weightedCol = weightedCol/totalWeight;
1477 weightedRow = weightedRow/totalWeight;
1478 geom->RecalculateTowerPosition(weightedRow, weightedCol, iSupModMax, depth, fMisalTransShift, fMisalRotShift, xyzNew);
1482 //printf("In Different SM\n");
1483 geom->RecalculateTowerPosition(iphi, ieta, iSupModMax, depth, fMisalTransShift, fMisalRotShift, xyzNew);
1486 clu->SetPosition(xyzNew);
1489 //___________________________________________________________________________________________
1490 void AliEMCALRecoUtils::RecalculateClusterDistanceToBadChannel(const AliEMCALGeometry * geom,
1491 AliVCaloCells* cells,
1492 AliVCluster * cluster)
1494 //re-evaluate distance to bad channel with updated bad map
1496 if(!fRecalDistToBadChannels) return;
1500 AliInfo("Cluster pointer null!");
1504 //Get channels map of the supermodule where the cluster is.
1505 Int_t absIdMax = -1, iSupMod =-1, icolM = -1, irowM = -1;
1506 Bool_t shared = kFALSE;
1507 GetMaxEnergyCell(geom, cells, cluster, absIdMax, iSupMod, icolM, irowM, shared);
1508 TH2D* hMap = (TH2D*)fEMCALBadChannelMap->At(iSupMod);
1511 Float_t minDist = 10000.;
1514 //Loop on tower status map
1515 for(Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++)
1517 for(Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++)
1519 //Check if tower is bad.
1520 if(hMap->GetBinContent(icol,irow)==0) continue;
1521 //printf("AliEMCALRecoUtils::RecalculateDistanceToBadChannels() - \n \t Bad channel in SM %d, col %d, row %d, \n \t Cluster max in col %d, row %d\n",
1522 // iSupMod,icol, irow, icolM,irowM);
1524 dRrow=TMath::Abs(irowM-irow);
1525 dRcol=TMath::Abs(icolM-icol);
1526 dist=TMath::Sqrt(dRrow*dRrow+dRcol*dRcol);
1529 //printf("MIN DISTANCE TO BAD %2.2f\n",dist);
1535 //In case the cluster is shared by 2 SuperModules, need to check the map of the second Super Module
1539 Int_t iSupMod2 = -1;
1541 //The only possible combinations are (0,1), (2,3) ... (8,9)
1542 if(iSupMod%2) iSupMod2 = iSupMod-1;
1543 else iSupMod2 = iSupMod+1;
1544 hMap2 = (TH2D*)fEMCALBadChannelMap->At(iSupMod2);
1546 //Loop on tower status map of second super module
1547 for(Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++)
1549 for(Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++)
1551 //Check if tower is bad.
1552 if(hMap2->GetBinContent(icol,irow)==0) continue;
1553 //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",
1554 // iSupMod2,icol, irow,iSupMod,icolM,irowM);
1555 dRrow=TMath::Abs(irow-irowM);
1559 dRcol=TMath::Abs(icol-(AliEMCALGeoParams::fgkEMCALCols+icolM));
1562 dRcol=TMath::Abs(AliEMCALGeoParams::fgkEMCALCols+icol-icolM);
1565 dist=TMath::Sqrt(dRrow*dRrow+dRcol*dRcol);
1566 if(dist < minDist) minDist = dist;
1569 }// shared cluster in 2 SuperModules
1571 AliDebug(2,Form("Max cluster cell (SM,col,row)=(%d %d %d) - Distance to Bad Channel %2.2f",iSupMod, icolM, irowM, minDist));
1572 cluster->SetDistanceToBadChannel(minDist);
1575 //__________________________________________________________________
1576 void AliEMCALRecoUtils::RecalculateClusterPID(AliVCluster * cluster)
1578 //re-evaluate identification parameters with bayesian
1582 AliInfo("Cluster pointer null!");
1586 if ( cluster->GetM02() != 0)
1587 fPIDUtils->ComputePID(cluster->E(),cluster->GetM02());
1589 Float_t pidlist[AliPID::kSPECIESCN+1];
1590 for(Int_t i = 0; i < AliPID::kSPECIESCN+1; i++) pidlist[i] = fPIDUtils->GetPIDFinal(i);
1592 cluster->SetPID(pidlist);
1595 //___________________________________________________________________________________________________________________
1596 void AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(const AliEMCALGeometry * geom,
1597 AliVCaloCells* cells,
1598 AliVCluster * cluster,
1599 Float_t & l0, Float_t & l1,
1600 Float_t & disp, Float_t & dEta, Float_t & dPhi,
1601 Float_t & sEta, Float_t & sPhi, Float_t & sEtaPhi)
1603 // Calculates new center of gravity in the local EMCAL-module coordinates
1604 // and tranfers into global ALICE coordinates
1605 // Calculates Dispersion and main axis
1609 AliInfo("Cluster pointer null!");
1613 Double_t eCell = 0.;
1614 Float_t fraction = 1.;
1615 Float_t recalFactor = 1.;
1623 Double_t etai = -1.;
1624 Double_t phii = -1.;
1629 Double_t etaMean = 0.;
1630 Double_t phiMean = 0.;
1632 //Loop on cells, calculate the cluster energy, in case a cut on cell energy is added
1633 // and to check if the cluster is between 2 SM in eta
1635 Bool_t shared = kFALSE;
1638 for(Int_t iDigit=0; iDigit < cluster->GetNCells(); iDigit++)
1640 //Get from the absid the supermodule, tower and eta/phi numbers
1641 geom->GetCellIndex(cluster->GetCellAbsId(iDigit),iSupMod,iTower,iIphi,iIeta);
1642 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi,iIeta, iphi,ieta);
1644 //Check if there are cells of different SM
1645 if (iDigit == 0 ) iSM0 = iSupMod;
1646 else if(iSupMod!= iSM0) shared = kTRUE;
1648 //Get the cell energy, if recalibration is on, apply factors
1649 fraction = cluster->GetCellAmplitudeFraction(iDigit);
1650 if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1652 if(IsRecalibrationOn())
1654 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1657 eCell = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1664 for(Int_t iDigit=0; iDigit < cluster->GetNCells(); iDigit++)
1666 //Get from the absid the supermodule, tower and eta/phi numbers
1667 geom->GetCellIndex(cluster->GetCellAbsId(iDigit),iSupMod,iTower,iIphi,iIeta);
1668 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi,iIeta, iphi,ieta);
1670 //Get the cell energy, if recalibration is on, apply factors
1671 fraction = cluster->GetCellAmplitudeFraction(iDigit);
1672 if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1674 if (!fCellsRecalibrated)
1676 if(IsRecalibrationOn())
1678 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1682 eCell = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1684 // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
1685 // C Side impair SM, nSupMod%2=1; A side pair SM, nSupMod%2=0
1686 if(shared && iSupMod%2) ieta+=AliEMCALGeoParams::fgkEMCALCols;
1688 if(cluster->E() > 0 && eCell > 0)
1690 w = GetCellWeight(eCell,cluster->E());
1692 etai=(Double_t)ieta;
1693 phii=(Double_t)iphi;
1700 sEta += w * etai * etai ;
1701 etaMean += w * etai ;
1702 sPhi += w * phii * phii ;
1703 phiMean += w * phii ;
1704 sEtaPhi += w * etai * phii ;
1708 AliError(Form("Wrong energy %f and/or amplitude %f\n", eCell, cluster->E()));
1711 //Normalize to the weight
1718 AliError(Form("Wrong weight %f\n", wtot));
1720 //Calculate dispersion
1721 for(Int_t iDigit=0; iDigit < cluster->GetNCells(); iDigit++)
1723 //Get from the absid the supermodule, tower and eta/phi numbers
1724 geom->GetCellIndex(cluster->GetCellAbsId(iDigit),iSupMod,iTower,iIphi,iIeta);
1725 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi,iIeta, iphi,ieta);
1727 //Get the cell energy, if recalibration is on, apply factors
1728 fraction = cluster->GetCellAmplitudeFraction(iDigit);
1729 if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1730 if (IsRecalibrationOn())
1732 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1734 eCell = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1736 // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
1737 // C Side impair SM, nSupMod%2=1; A side pair SM, nSupMod%2=0
1738 if(shared && iSupMod%2) ieta+=AliEMCALGeoParams::fgkEMCALCols;
1740 if(cluster->E() > 0 && eCell > 0)
1742 w = GetCellWeight(eCell,cluster->E());
1744 etai=(Double_t)ieta;
1745 phii=(Double_t)iphi;
1748 disp += w *((etai-etaMean)*(etai-etaMean)+(phii-phiMean)*(phii-phiMean));
1749 dEta += w * (etai-etaMean)*(etai-etaMean) ;
1750 dPhi += w * (phii-phiMean)*(phii-phiMean) ;
1754 AliError(Form("Wrong energy %f and/or amplitude %f\n", eCell, cluster->E()));
1757 //Normalize to the weigth and set shower shape parameters
1758 if (wtot > 0 && nstat > 1)
1767 sEta -= etaMean * etaMean ;
1768 sPhi -= phiMean * phiMean ;
1769 sEtaPhi -= etaMean * phiMean ;
1771 l0 = (0.5 * (sEta + sPhi) + TMath::Sqrt( 0.25 * (sEta - sPhi) * (sEta - sPhi) + sEtaPhi * sEtaPhi ));
1772 l1 = (0.5 * (sEta + sPhi) - TMath::Sqrt( 0.25 * (sEta - sPhi) * (sEta - sPhi) + sEtaPhi * sEtaPhi ));
1778 dEta = 0. ; dPhi = 0. ; disp = 0. ;
1779 sEta = 0. ; sPhi = 0. ; sEtaPhi = 0. ;
1784 //____________________________________________________________________________________________
1785 void AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(const AliEMCALGeometry * geom,
1786 AliVCaloCells* cells,
1787 AliVCluster * cluster)
1789 // Calculates new center of gravity in the local EMCAL-module coordinates
1790 // and tranfers into global ALICE coordinates
1791 // Calculates Dispersion and main axis and puts them into the cluster
1793 Float_t l0 = 0., l1 = 0.;
1794 Float_t disp = 0., dEta = 0., dPhi = 0.;
1795 Float_t sEta = 0., sPhi = 0., sEtaPhi = 0.;
1797 AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(geom,cells,cluster,l0,l1,disp,
1798 dEta, dPhi, sEta, sPhi, sEtaPhi);
1800 cluster->SetM02(l0);
1801 cluster->SetM20(l1);
1802 if(disp > 0. ) cluster->SetDispersion(TMath::Sqrt(disp)) ;
1806 //____________________________________________________________________________
1807 void AliEMCALRecoUtils::FindMatches(AliVEvent *event,
1808 TObjArray * clusterArr,
1809 const AliEMCALGeometry *geom)
1811 //This function should be called before the cluster loop
1812 //Before call this function, please recalculate the cluster positions
1813 //Given the input event, loop over all the tracks, select the closest cluster as matched with fCutR
1814 //Store matched cluster indexes and residuals
1816 fMatchedTrackIndex ->Reset();
1817 fMatchedClusterIndex->Reset();
1818 fResidualPhi->Reset();
1819 fResidualEta->Reset();
1821 fMatchedTrackIndex ->Set(1000);
1822 fMatchedClusterIndex->Set(1000);
1823 fResidualPhi->Set(1000);
1824 fResidualEta->Set(1000);
1826 AliESDEvent* esdevent = dynamic_cast<AliESDEvent*> (event);
1827 AliAODEvent* aodevent = dynamic_cast<AliAODEvent*> (event);
1829 // init the magnetic field if not already on
1830 if(!TGeoGlobalMagField::Instance()->GetField())
1832 if (!event->InitMagneticField())
1834 AliInfo("Mag Field not initialized, null esd/aod evetn pointers");
1840 UInt_t mask1 = esdevent->GetESDRun()->GetDetectorsInDAQ();
1841 UInt_t mask2 = esdevent->GetESDRun()->GetDetectorsInReco();
1842 Bool_t desc1 = (mask1 >> 3) & 0x1;
1843 Bool_t desc2 = (mask2 >> 3) & 0x1;
1844 if (desc1==0 || desc2==0) {
1845 // AliError(Form("TPC not in DAQ/RECO: %u (%u)/%u (%u)",
1846 // mask1, esdevent->GetESDRun()->GetDetectorsInReco(),
1847 // mask2, esdevent->GetESDRun()->GetDetectorsInDAQ()));
1852 TObjArray *clusterArray = 0x0;
1855 clusterArray = new TObjArray(event->GetNumberOfCaloClusters());
1856 for(Int_t icl=0; icl<event->GetNumberOfCaloClusters(); icl++)
1858 AliVCluster *cluster = (AliVCluster*) event->GetCaloCluster(icl);
1859 if(geom && !IsGoodCluster(cluster,geom,(AliVCaloCells*)event->GetEMCALCells())) continue;
1860 clusterArray->AddAt(cluster,icl);
1866 for (Int_t i=0; i<21;i++) cv[i]=0;
1867 for(Int_t itr=0; itr<event->GetNumberOfTracks(); itr++)
1869 AliExternalTrackParam *trackParam = 0;
1871 //If the input event is ESD, the starting point for extrapolation is TPCOut, if available, or TPCInner
1872 AliESDtrack *esdTrack = 0;
1873 AliAODTrack *aodTrack = 0;
1876 esdTrack = esdevent->GetTrack(itr);
1877 if(!esdTrack) continue;
1878 if(!IsAccepted(esdTrack)) continue;
1879 if(esdTrack->Pt()<fCutMinTrackPt) continue;
1880 Double_t phi = esdTrack->Phi()*TMath::RadToDeg();
1881 if(TMath::Abs(esdTrack->Eta())>0.8 || phi <= 20 || phi >= 240 ) continue;
1883 trackParam = const_cast<AliExternalTrackParam*>(esdTrack->GetInnerParam()); // if TPC Available
1885 trackParam = new AliExternalTrackParam(*esdTrack); // If ITS Track Standing alone
1888 //If the input event is AOD, the starting point for extrapolation is at vertex
1889 //AOD tracks are selected according to its filterbit.
1892 aodTrack = aodevent->GetTrack(itr);
1893 if(!aodTrack) continue;
1895 if(fAODTPCOnlyTracks) // Match with TPC only tracks, default from May 2013, before filter bit 32
1897 //printf("Match with TPC only tracks, accept? %d, test bit 128 <%d> \n", aodTrack->IsTPCOnly(), aodTrack->TestFilterMask(128));
1898 if(!aodTrack->IsTPCOnly()) continue ;
1900 else if(fAODHybridTracks) // Match with hybrid tracks
1902 //printf("Match with Hybrid tracks, accept? %d \n", aodTrack->IsHybridGlobalConstrainedGlobal());
1903 if(!aodTrack->IsHybridGlobalConstrainedGlobal()) continue ;
1905 else // Match with tracks on a mask
1907 //printf("Match with tracks having filter bit mask %d, accept? %d \n",fAODFilterMask,aodTrack->TestFilterMask(fAODFilterMask));
1908 if(!aodTrack->TestFilterMask(fAODFilterMask) ) continue; //Select AOD tracks
1911 if(aodTrack->Pt()<fCutMinTrackPt) continue;
1913 Double_t phi = aodTrack->Phi()*TMath::RadToDeg();
1914 if(TMath::Abs(aodTrack->Eta())>0.8 || phi <= 20 || phi >= 240 ) continue;
1915 Double_t pos[3],mom[3];
1916 aodTrack->GetXYZ(pos);
1917 aodTrack->GetPxPyPz(mom);
1918 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()));
1920 trackParam= new AliExternalTrackParam(pos,mom,cv,aodTrack->Charge());
1923 //Return if the input data is not "AOD" or "ESD"
1926 printf("Wrong input data type! Should be \"AOD\" or \"ESD\"\n");
1929 clusterArray->Clear();
1930 delete clusterArray;
1935 if(!trackParam) continue;
1937 //Extrapolate the track to EMCal surface
1938 AliExternalTrackParam emcalParam(*trackParam);
1939 Float_t eta, phi, pt;
1940 if(!ExtrapolateTrackToEMCalSurface(&emcalParam, fEMCalSurfaceDistance, fMass, fStepSurface, eta, phi, pt))
1942 if(aodevent && trackParam) delete trackParam;
1943 if(fITSTrackSA && trackParam) delete trackParam;
1947 if(TMath::Abs(eta)>0.75 || (phi) < 70*TMath::DegToRad() || (phi) > 190*TMath::DegToRad())
1949 if(aodevent && trackParam) delete trackParam;
1950 if(fITSTrackSA && trackParam) delete trackParam;
1954 //Find matched clusters
1956 Float_t dEta = -999, dPhi = -999;
1959 index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArray, dEta, dPhi);
1963 index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArr, dEta, dPhi);
1968 fMatchedTrackIndex ->AddAt(itr,matched);
1969 fMatchedClusterIndex ->AddAt(index,matched);
1970 fResidualEta ->AddAt(dEta,matched);
1971 fResidualPhi ->AddAt(dPhi,matched);
1974 if(aodevent && trackParam) delete trackParam;
1975 if(fITSTrackSA && trackParam) delete trackParam;
1980 clusterArray->Clear();
1981 delete clusterArray;
1984 AliDebug(2,Form("Number of matched pairs = %d !\n",matched));
1986 fMatchedTrackIndex ->Set(matched);
1987 fMatchedClusterIndex ->Set(matched);
1988 fResidualPhi ->Set(matched);
1989 fResidualEta ->Set(matched);
1992 //________________________________________________________________________________
1993 Int_t AliEMCALRecoUtils::FindMatchedClusterInEvent(const AliESDtrack *track,
1994 const AliVEvent *event,
1995 const AliEMCALGeometry *geom,
1996 Float_t &dEta, Float_t &dPhi)
1999 // This function returns the index of matched cluster to input track
2000 // Returns -1 if no match is found
2002 Double_t phiV = track->Phi()*TMath::RadToDeg();
2003 if(TMath::Abs(track->Eta())>0.8 || phiV <= 20 || phiV >= 240 ) return index;
2004 AliExternalTrackParam *trackParam = 0;
2006 trackParam = const_cast<AliExternalTrackParam*>(track->GetInnerParam()); // If TPC
2008 trackParam = new AliExternalTrackParam(*track);
2010 if(!trackParam) return index;
2011 AliExternalTrackParam emcalParam(*trackParam);
2012 Float_t eta, phi, pt;
2014 if(!ExtrapolateTrackToEMCalSurface(&emcalParam, fEMCalSurfaceDistance, fMass, fStepSurface, eta, phi, pt)) {
2015 if(fITSTrackSA) delete trackParam;
2018 if(TMath::Abs(eta)>0.75 || (phi) < 70*TMath::DegToRad() || (phi) > 190*TMath::DegToRad()){
2019 if(fITSTrackSA) delete trackParam;
2023 TObjArray *clusterArr = new TObjArray(event->GetNumberOfCaloClusters());
2025 for(Int_t icl=0; icl<event->GetNumberOfCaloClusters(); icl++)
2027 AliVCluster *cluster = (AliVCluster*) event->GetCaloCluster(icl);
2028 if(geom && !IsGoodCluster(cluster,geom,(AliVCaloCells*)event->GetEMCALCells())) continue;
2029 clusterArr->AddAt(cluster,icl);
2032 index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArr, dEta, dPhi);
2033 clusterArr->Clear();
2035 if(fITSTrackSA) delete trackParam;
2040 //_______________________________________________________________________________________________
2041 Int_t AliEMCALRecoUtils::FindMatchedClusterInClusterArr(const AliExternalTrackParam *emcalParam,
2042 AliExternalTrackParam *trkParam,
2043 const TObjArray * clusterArr,
2044 Float_t &dEta, Float_t &dPhi)
2046 // Find matched cluster in array
2048 dEta=-999, dPhi=-999;
2049 Float_t dRMax = fCutR, dEtaMax=fCutEta, dPhiMax=fCutPhi;
2051 Float_t tmpEta=-999, tmpPhi=-999;
2053 Double_t exPos[3] = {0.,0.,0.};
2054 if(!emcalParam->GetXYZ(exPos)) return index;
2056 Float_t clsPos[3] = {0.,0.,0.};
2057 for(Int_t icl=0; icl<clusterArr->GetEntriesFast(); icl++)
2059 AliVCluster *cluster = dynamic_cast<AliVCluster*> (clusterArr->At(icl)) ;
2060 if(!cluster || !cluster->IsEMCAL()) continue;
2061 cluster->GetPosition(clsPos);
2062 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));
2063 if(dR > fClusterWindow) continue;
2065 AliExternalTrackParam trkPamTmp (*trkParam);//Retrieve the starting point every time before the extrapolation
2066 if(!ExtrapolateTrackToCluster(&trkPamTmp, cluster, fMass, fStepCluster, tmpEta, tmpPhi)) continue;
2069 Float_t tmpR=TMath::Sqrt(tmpEta*tmpEta + tmpPhi*tmpPhi);
2078 else if(fCutEtaPhiSeparate)
2080 if(TMath::Abs(tmpEta)<TMath::Abs(dEtaMax) && TMath::Abs(tmpPhi)<TMath::Abs(dPhiMax))
2089 printf("Error: please specify your cut criteria\n");
2090 printf("To cut on sqrt(dEta^2+dPhi^2), use: SwitchOnCutEtaPhiSum()\n");
2091 printf("To cut on dEta and dPhi separately, use: SwitchOnCutEtaPhiSeparate()\n");
2102 //------------------------------------------------------------------------------------
2103 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToEMCalSurface(AliVTrack *track,
2104 Double_t emcalR, Double_t mass, Double_t step)
2106 // Extrpolate track to EMCAL surface
2108 track->SetTrackPhiEtaPtOnEMCal(-999, -999, -999);
2110 if (track->Pt()<0.350)
2113 Double_t phi = track->Phi()*TMath::RadToDeg();
2114 if (TMath::Abs(track->Eta())>0.9 || phi <= 10 || phi >= 250)
2117 AliESDtrack *esdt = dynamic_cast<AliESDtrack*>(track);
2118 AliAODTrack *aodt = 0;
2120 aodt = dynamic_cast<AliAODTrack*>(track);
2126 Bool_t onlyTPC = kFALSE;
2130 mass = esdt->GetMass(onlyTPC);
2135 AliExternalTrackParam *trackParam = 0;
2137 trackParam = new AliExternalTrackParam(*esdt->GetInnerParam());
2139 Double_t xyz[3] = {0}, pxpypz[3] = {0}, cv[21] = {0};
2140 aodt->PxPyPz(pxpypz);
2142 aodt->GetCovarianceXYZPxPyPz(cv);
2143 trackParam = new AliExternalTrackParam(xyz,pxpypz,cv,aodt->Charge());
2148 Float_t etaout=-999, phiout=-999, ptout=-999;
2149 Bool_t ret = ExtrapolateTrackToEMCalSurface(trackParam,
2159 if (TMath::Abs(etaout)>0.75 || (phiout<70*TMath::DegToRad()) || (phiout>190*TMath::DegToRad()))
2161 track->SetTrackPhiEtaPtOnEMCal(phiout, etaout, ptout);
2166 //------------------------------------------------------------------------------------
2167 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToEMCalSurface(AliExternalTrackParam *trkParam,
2168 const Double_t emcalR,
2169 const Double_t mass,
2170 const Double_t step,
2175 //Extrapolate track to EMCAL surface
2177 eta = -999, phi = -999, pt = -999;
2178 if(!trkParam) return kFALSE;
2179 if(!AliTrackerBase::PropagateTrackToBxByBz(trkParam, emcalR, mass, step, kTRUE, 0.8, -1)) return kFALSE;
2180 Double_t trkPos[3] = {0.,0.,0.};
2181 if(!trkParam->GetXYZ(trkPos)) return kFALSE;
2182 TVector3 trkPosVec(trkPos[0],trkPos[1],trkPos[2]);
2183 eta = trkPosVec.Eta();
2184 phi = trkPosVec.Phi();
2185 pt = trkParam->Pt();
2187 phi += 2*TMath::Pi();
2192 //-----------------------------------------------------------------------------------
2193 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToPosition(AliExternalTrackParam *trkParam,
2194 const Float_t *clsPos,
2201 //Return the residual by extrapolating a track param to a global position
2205 if(!trkParam) return kFALSE;
2206 Double_t trkPos[3] = {0.,0.,0.};
2207 TVector3 vec(clsPos[0],clsPos[1],clsPos[2]);
2208 Double_t alpha = ((int)(vec.Phi()*TMath::RadToDeg()/20)+0.5)*20*TMath::DegToRad();
2209 vec.RotateZ(-alpha); //Rotate the cluster to the local extrapolation coordinate system
2210 if(!AliTrackerBase::PropagateTrackToBxByBz(trkParam, vec.X(), mass, step,kTRUE, 0.8, -1)) return kFALSE;
2211 if(!trkParam->GetXYZ(trkPos)) return kFALSE; //Get the extrapolated global position
2213 TVector3 clsPosVec(clsPos[0],clsPos[1],clsPos[2]);
2214 TVector3 trkPosVec(trkPos[0],trkPos[1],trkPos[2]);
2216 // track cluster matching
2217 tmpPhi = clsPosVec.DeltaPhi(trkPosVec); // tmpPhi is between -pi and pi
2218 tmpEta = clsPosVec.Eta()-trkPosVec.Eta();
2223 //----------------------------------------------------------------------------------
2224 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToCluster(AliExternalTrackParam *trkParam,
2225 const AliVCluster *cluster,
2226 const Double_t mass,
2227 const Double_t step,
2232 //Return the residual by extrapolating a track param to a cluster
2236 if(!cluster || !trkParam) return kFALSE;
2238 Float_t clsPos[3] = {0.,0.,0.};
2239 cluster->GetPosition(clsPos);
2241 return ExtrapolateTrackToPosition(trkParam, clsPos, mass, step, tmpEta, tmpPhi);
2244 //---------------------------------------------------------------------------------
2245 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToCluster(AliExternalTrackParam *trkParam,
2246 const AliVCluster *cluster,
2251 //Return the residual by extrapolating a track param to a clusterfStepCluster
2254 return ExtrapolateTrackToCluster(trkParam, cluster, fMass, fStepCluster, tmpEta, tmpPhi);
2257 //_______________________________________________________________________
2258 void AliEMCALRecoUtils::GetMatchedResiduals(const Int_t clsIndex,
2259 Float_t &dEta, Float_t &dPhi)
2261 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2262 //Get the residuals dEta and dPhi for this cluster to the closest track
2263 //Works with ESDs and AODs
2265 if( FindMatchedPosForCluster(clsIndex) >= 999 )
2267 AliDebug(2,"No matched tracks found!\n");
2272 dEta = fResidualEta->At(FindMatchedPosForCluster(clsIndex));
2273 dPhi = fResidualPhi->At(FindMatchedPosForCluster(clsIndex));
2276 //______________________________________________________________________________________________
2277 void AliEMCALRecoUtils::GetMatchedClusterResiduals(Int_t trkIndex, Float_t &dEta, Float_t &dPhi)
2279 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2280 //Get the residuals dEta and dPhi for this track to the closest cluster
2281 //Works with ESDs and AODs
2283 if( FindMatchedPosForTrack(trkIndex) >= 999 )
2285 AliDebug(2,"No matched cluster found!\n");
2290 dEta = fResidualEta->At(FindMatchedPosForTrack(trkIndex));
2291 dPhi = fResidualPhi->At(FindMatchedPosForTrack(trkIndex));
2294 //__________________________________________________________
2295 Int_t AliEMCALRecoUtils::GetMatchedTrackIndex(Int_t clsIndex)
2297 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2298 //Get the index of matched track to this cluster
2299 //Works with ESDs and AODs
2301 if(IsClusterMatched(clsIndex))
2302 return fMatchedTrackIndex->At(FindMatchedPosForCluster(clsIndex));
2307 //__________________________________________________________
2308 Int_t AliEMCALRecoUtils::GetMatchedClusterIndex(Int_t trkIndex)
2310 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2311 //Get the index of matched cluster to this track
2312 //Works with ESDs and AODs
2314 if(IsTrackMatched(trkIndex))
2315 return fMatchedClusterIndex->At(FindMatchedPosForTrack(trkIndex));
2320 //______________________________________________________________
2321 Bool_t AliEMCALRecoUtils::IsClusterMatched(Int_t clsIndex) const
2323 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2324 //Returns if the cluster has a match
2325 if(FindMatchedPosForCluster(clsIndex) < 999)
2331 //____________________________________________________________
2332 Bool_t AliEMCALRecoUtils::IsTrackMatched(Int_t trkIndex) const
2334 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2335 //Returns if the track has a match
2336 if(FindMatchedPosForTrack(trkIndex) < 999)
2342 //______________________________________________________________________
2343 UInt_t AliEMCALRecoUtils::FindMatchedPosForCluster(Int_t clsIndex) const
2345 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2346 //Returns the position of the match in the fMatchedClusterIndex array
2347 Float_t tmpR = fCutR;
2350 for(Int_t i=0; i<fMatchedClusterIndex->GetSize(); i++)
2352 if(fMatchedClusterIndex->At(i)==clsIndex)
2354 Float_t r = TMath::Sqrt(fResidualEta->At(i)*fResidualEta->At(i) + fResidualPhi->At(i)*fResidualPhi->At(i));
2359 AliDebug(3,Form("Matched cluster index: index: %d, dEta: %2.4f, dPhi: %2.4f.\n",
2360 fMatchedClusterIndex->At(i),fResidualEta->At(i),fResidualPhi->At(i)));
2367 //____________________________________________________________________
2368 UInt_t AliEMCALRecoUtils::FindMatchedPosForTrack(Int_t trkIndex) const
2370 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2371 //Returns the position of the match in the fMatchedTrackIndex array
2372 Float_t tmpR = fCutR;
2375 for(Int_t i=0; i<fMatchedTrackIndex->GetSize(); i++)
2377 if(fMatchedTrackIndex->At(i)==trkIndex)
2379 Float_t r = TMath::Sqrt(fResidualEta->At(i)*fResidualEta->At(i) + fResidualPhi->At(i)*fResidualPhi->At(i));
2384 AliDebug(3,Form("Matched track index: index: %d, dEta: %2.4f, dPhi: %2.4f.\n",
2385 fMatchedTrackIndex->At(i),fResidualEta->At(i),fResidualPhi->At(i)));
2392 //__________________________________________________________________________
2393 Bool_t AliEMCALRecoUtils::IsGoodCluster(AliVCluster *cluster,
2394 const AliEMCALGeometry *geom,
2395 AliVCaloCells* cells,const Int_t bc)
2397 // check if the cluster survives some quality cut
2400 Bool_t isGood=kTRUE;
2402 if(!cluster || !cluster->IsEMCAL()) return kFALSE;
2404 if(ClusterContainsBadChannel(geom,cluster->GetCellsAbsId(),cluster->GetNCells())) return kFALSE;
2406 if(!CheckCellFiducialRegion(geom,cluster,cells)) return kFALSE;
2408 if(IsExoticCluster(cluster, cells,bc)) return kFALSE;
2413 //__________________________________________________________
2414 Bool_t AliEMCALRecoUtils::IsAccepted(AliESDtrack *esdTrack)
2416 // Given a esd track, return whether the track survive all the cuts
2418 // The different quality parameter are first
2419 // retrieved from the track. then it is found out what cuts the
2420 // track did not survive and finally the cuts are imposed.
2422 UInt_t status = esdTrack->GetStatus();
2424 Int_t nClustersITS = esdTrack->GetITSclusters(0);
2425 Int_t nClustersTPC = esdTrack->GetTPCclusters(0);
2427 Float_t chi2PerClusterITS = -1;
2428 Float_t chi2PerClusterTPC = -1;
2429 if (nClustersITS!=0)
2430 chi2PerClusterITS = esdTrack->GetITSchi2()/Float_t(nClustersITS);
2431 if (nClustersTPC!=0)
2432 chi2PerClusterTPC = esdTrack->GetTPCchi2()/Float_t(nClustersTPC);
2436 if(fTrackCutsType==kGlobalCut)
2438 Float_t maxDCAToVertexXYPtDep = 0.0182 + 0.0350/TMath::Power(esdTrack->Pt(),1.01); //This expression comes from AliESDtrackCuts::GetStandardITSTPCTrackCuts2010()
2439 //AliDebug(3,Form("Track pT = %f, DCAtoVertexXY = %f",esdTrack->Pt(),MaxDCAToVertexXYPtDep));
2440 SetMaxDCAToVertexXY(maxDCAToVertexXYPtDep); //Set pT dependent DCA cut to vertex in x-y plane
2446 esdTrack->GetImpactParameters(b,bCov);
2447 if (bCov[0]<=0 || bCov[2]<=0)
2449 AliDebug(1, "Estimated b resolution lower or equal zero!");
2450 bCov[0]=0; bCov[2]=0;
2453 Float_t dcaToVertexXY = b[0];
2454 Float_t dcaToVertexZ = b[1];
2455 Float_t dcaToVertex = -1;
2457 if (fCutDCAToVertex2D)
2458 dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY/fCutMaxDCAToVertexXY/fCutMaxDCAToVertexXY + dcaToVertexZ*dcaToVertexZ/fCutMaxDCAToVertexZ/fCutMaxDCAToVertexZ);
2460 dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY + dcaToVertexZ*dcaToVertexZ);
2464 Bool_t cuts[kNCuts];
2465 for (Int_t i=0; i<kNCuts; i++) cuts[i]=kFALSE;
2467 // track quality cuts
2468 if (fCutRequireTPCRefit && (status&AliESDtrack::kTPCrefit)==0)
2470 if (fCutRequireITSRefit && (status&AliESDtrack::kITSrefit)==0)
2472 if (nClustersTPC<fCutMinNClusterTPC)
2474 if (nClustersITS<fCutMinNClusterITS)
2476 if (chi2PerClusterTPC>fCutMaxChi2PerClusterTPC)
2478 if (chi2PerClusterITS>fCutMaxChi2PerClusterITS)
2480 if (!fCutAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0)
2482 if (fCutDCAToVertex2D && dcaToVertex > 1)
2484 if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexXY) > fCutMaxDCAToVertexXY)
2486 if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexZ) > fCutMaxDCAToVertexZ)
2489 if(fTrackCutsType==kGlobalCut)
2491 //Require at least one SPD point + anything else in ITS
2492 if( (esdTrack->HasPointOnITSLayer(0) || esdTrack->HasPointOnITSLayer(1)) == kFALSE)
2497 if(fCutRequireITSStandAlone || fCutRequireITSpureSA){
2498 if ((status & AliESDtrack::kITSin) == 0 || (status & AliESDtrack::kTPCin)){
2502 // ITS standalone tracks
2503 if(fCutRequireITSStandAlone && !fCutRequireITSpureSA){
2504 if(status & AliESDtrack::kITSpureSA) cuts[11] = kTRUE;
2505 }else if(fCutRequireITSpureSA){
2506 if(!(status & AliESDtrack::kITSpureSA)) cuts[11] = kTRUE;
2512 for (Int_t i=0; i<kNCuts; i++)
2513 if (cuts[i]) { cut = kTRUE ; }
2522 //_____________________________________
2523 void AliEMCALRecoUtils::InitTrackCuts()
2525 //Intilize the track cut criteria
2526 //By default these cuts are set according to AliESDtrackCuts::GetStandardTPCOnlyTrackCuts()
2527 //Also you can customize the cuts using the setters
2529 switch (fTrackCutsType)
2533 AliInfo(Form("Track cuts for matching: GetStandardTPCOnlyTrackCuts()"));
2535 SetMinNClustersTPC(70);
2536 SetMaxChi2PerClusterTPC(4);
2537 SetAcceptKinkDaughters(kFALSE);
2538 SetRequireTPCRefit(kFALSE);
2541 SetRequireITSRefit(kFALSE);
2542 SetMaxDCAToVertexZ(3.2);
2543 SetMaxDCAToVertexXY(2.4);
2544 SetDCAToVertex2D(kTRUE);
2551 AliInfo(Form("Track cuts for matching: GetStandardITSTPCTrackCuts2010(kTURE)"));
2553 SetMinNClustersTPC(70);
2554 SetMaxChi2PerClusterTPC(4);
2555 SetAcceptKinkDaughters(kFALSE);
2556 SetRequireTPCRefit(kTRUE);
2559 SetRequireITSRefit(kTRUE);
2560 SetMaxDCAToVertexZ(2);
2561 SetMaxDCAToVertexXY();
2562 SetDCAToVertex2D(kFALSE);
2569 AliInfo(Form("Track cuts for matching: Loose cut w/o DCA cut"));
2570 SetMinNClustersTPC(50);
2571 SetAcceptKinkDaughters(kTRUE);
2576 case kITSStandAlone:
2578 AliInfo(Form("Track cuts for matching: ITS Stand Alone tracks cut w/o DCA cut"));
2579 SetRequireITSRefit(kTRUE);
2580 SetRequireITSStandAlone(kTRUE);
2581 SetITSTrackSA(kTRUE);
2589 //________________________________________________________________________
2590 void AliEMCALRecoUtils::SetClusterMatchedToTrack(const AliVEvent *event)
2592 // Checks if tracks are matched to EMC clusters and set the matched EMCAL cluster index to ESD track.
2594 Int_t nTracks = event->GetNumberOfTracks();
2595 for (Int_t iTrack = 0; iTrack < nTracks; ++iTrack)
2597 AliVTrack* track = dynamic_cast<AliVTrack*>(event->GetTrack(iTrack));
2600 AliWarning(Form("Could not receive track %d", iTrack));
2604 Int_t matchClusIndex = GetMatchedClusterIndex(iTrack);
2605 track->SetEMCALcluster(matchClusIndex); //sets -1 if track not matched within residual
2606 /*the following can be done better if AliVTrack::SetStatus will be there. Patch pending with Andreas/Peter*/
2607 AliESDtrack* esdtrack = dynamic_cast<AliESDtrack*>(track);
2609 if(matchClusIndex != -1)
2610 esdtrack->SetStatus(AliESDtrack::kEMCALmatch);
2612 esdtrack->ResetStatus(AliESDtrack::kEMCALmatch);
2614 AliAODTrack* aodtrack = dynamic_cast<AliAODTrack*>(track);
2615 if(matchClusIndex != -1)
2616 aodtrack->SetStatus(AliESDtrack::kEMCALmatch);
2618 aodtrack->ResetStatus(AliESDtrack::kEMCALmatch);
2622 AliDebug(2,"Track matched to closest cluster");
2625 //_________________________________________________________________________
2626 void AliEMCALRecoUtils::SetTracksMatchedToCluster(const AliVEvent *event)
2628 // Checks if EMC clusters are matched to ESD track.
2629 // Adds track indexes of all the tracks matched to a cluster withing residuals in ESDCalocluster.
2631 for (Int_t iClus=0; iClus < event->GetNumberOfCaloClusters(); ++iClus)
2633 AliVCluster *cluster = event->GetCaloCluster(iClus);
2634 if (!cluster->IsEMCAL())
2637 Int_t nTracks = event->GetNumberOfTracks();
2638 TArrayI arrayTrackMatched(nTracks);
2640 // Get the closest track matched to the cluster
2642 Int_t matchTrackIndex = GetMatchedTrackIndex(iClus);
2643 if (matchTrackIndex != -1)
2645 arrayTrackMatched[nMatched] = matchTrackIndex;
2649 // Get all other tracks matched to the cluster
2650 for(Int_t iTrk=0; iTrk<nTracks; ++iTrk)
2652 AliVTrack* track = dynamic_cast<AliVTrack*>(event->GetTrack(iTrk));
2653 if(iTrk == matchTrackIndex) continue;
2654 if(track->GetEMCALcluster() == iClus)
2656 arrayTrackMatched[nMatched] = iTrk;
2661 //printf("Tender::SetTracksMatchedToCluster - cluster E %f, N matches %d, first match %d\n",cluster->E(),nMatched,arrayTrackMatched[0]);
2663 arrayTrackMatched.Set(nMatched);
2664 AliESDCaloCluster *esdcluster = dynamic_cast<AliESDCaloCluster*>(cluster);
2666 esdcluster->AddTracksMatched(arrayTrackMatched);
2667 else if (nMatched>0) {
2668 AliAODCaloCluster *aodcluster = dynamic_cast<AliAODCaloCluster*>(cluster);
2670 aodcluster->AddTrackMatched(event->GetTrack(arrayTrackMatched.At(0)));
2673 Float_t eta= -999, phi = -999;
2674 if (matchTrackIndex != -1)
2675 GetMatchedResiduals(iClus, eta, phi);
2676 cluster->SetTrackDistance(phi, eta);
2679 AliDebug(2,"Cluster matched to tracks");
2682 //___________________________________________________
2683 void AliEMCALRecoUtils::Print(const Option_t *) const
2687 printf("AliEMCALRecoUtils Settings: \n");
2688 printf("Misalignment shifts\n");
2689 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,
2690 fMisalTransShift[i*3],fMisalTransShift[i*3+1],fMisalTransShift[i*3+2],
2691 fMisalRotShift[i*3], fMisalRotShift[i*3+1], fMisalRotShift[i*3+2] );
2692 printf("Non linearity function %d, parameters:\n", fNonLinearityFunction);
2693 for(Int_t i=0; i<6; i++) printf("param[%d]=%f\n",i, fNonLinearityParams[i]);
2695 printf("Position Recalculation option %d, Particle Type %d, fW0 %2.2f, Recalibrate Data %d \n",fPosAlgo,fParticleType,fW0, fRecalibration);
2697 printf("Matching criteria: ");
2700 printf("sqrt(dEta^2+dPhi^2)<%4.3f\n",fCutR);
2702 else if(fCutEtaPhiSeparate)
2704 printf("dEta<%4.3f, dPhi<%4.3f\n",fCutEta,fCutPhi);
2709 printf("please specify your cut criteria\n");
2710 printf("To cut on sqrt(dEta^2+dPhi^2), use: SwitchOnCutEtaPhiSum()\n");
2711 printf("To cut on dEta and dPhi separately, use: SwitchOnCutEtaPhiSeparate()\n");
2714 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);
2715 printf("Cluster selection window: dR < %2.0f\n",fClusterWindow);
2717 printf("Track cuts: \n");
2718 printf("Minimum track pT: %1.2f\n",fCutMinTrackPt);
2719 printf("AOD track selection: tpc only %d, or hybrid %d, or mask: %d\n",fAODTPCOnlyTracks,fAODHybridTracks, fAODFilterMask);
2720 printf("TPCRefit = %d, ITSRefit = %d\n",fCutRequireTPCRefit,fCutRequireITSRefit);
2721 printf("AcceptKinks = %d\n",fCutAcceptKinkDaughters);
2722 printf("MinNCulsterTPC = %d, MinNClusterITS = %d\n",fCutMinNClusterTPC,fCutMinNClusterITS);
2723 printf("MaxChi2TPC = %2.2f, MaxChi2ITS = %2.2f\n",fCutMaxChi2PerClusterTPC,fCutMaxChi2PerClusterITS);
2724 printf("DCSToVertex2D = %d, MaxDCAToVertexXY = %2.2f, MaxDCAToVertexZ = %2.2f\n",fCutDCAToVertex2D,fCutMaxDCAToVertexXY,fCutMaxDCAToVertexZ);