1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
16 /* $Id: AliEMCALRecoUtils.cxx 33808 2009-07-15 09:48:08Z gconesab $ */
18 ///////////////////////////////////////////////////////////////////////////////
20 // Class AliEMCALRecoUtils
21 // Some utilities to recalculate the cluster position or energy linearity
24 // Author: Gustavo Conesa (LPSC- Grenoble)
25 // Track matching part: Rongrong Ma (Yale)
27 ///////////////////////////////////////////////////////////////////////////////
30 // standard C++ includes
31 //#include <Riostream.h>
34 #include <TGeoManager.h>
35 #include <TGeoMatrix.h>
40 #include <TObjArray.h>
43 #include "AliVCluster.h"
44 #include "AliVCaloCells.h"
47 #include "AliESDEvent.h"
48 #include "AliAODEvent.h"
49 #include "AliESDtrack.h"
50 #include "AliAODTrack.h"
51 #include "AliExternalTrackParam.h"
52 #include "AliESDfriendTrack.h"
53 #include "AliTrackerBase.h"
56 #include "AliEMCALRecoUtils.h"
57 #include "AliEMCALGeometry.h"
58 #include "AliTrackerBase.h"
59 #include "AliEMCALCalibTimeDepCorrection.h" // Run dependent
60 #include "AliEMCALPIDUtils.h"
63 ClassImp(AliEMCALRecoUtils)
65 //_____________________________________
66 AliEMCALRecoUtils::AliEMCALRecoUtils():
67 fParticleType(0), fPosAlgo(0), fW0(0),
68 fNonLinearityFunction(0), fNonLinearThreshold(0),
69 fSmearClusterEnergy(kFALSE), fRandom(),
70 fCellsRecalibrated(kFALSE), fRecalibration(kFALSE), fEMCALRecalibrationFactors(),
71 fTimeRecalibration(kFALSE), fEMCALTimeRecalibrationFactors(),
72 fUseRunCorrectionFactors(kFALSE), fRunCorrectionFactorsSet(kFALSE),
73 fRemoveBadChannels(kFALSE), fRecalDistToBadChannels(kFALSE), fEMCALBadChannelMap(),
74 fNCellsFromEMCALBorder(0), fNoEMCALBorderAtEta0(kTRUE),
75 fRejectExoticCluster(kFALSE), fRejectExoticCells(kFALSE),
76 fExoticCellFraction(0), fExoticCellDiffTime(0), fExoticCellMinAmplitude(0),
77 fPIDUtils(), fAODFilterMask(0),
78 fMatchedTrackIndex(0x0), fMatchedClusterIndex(0x0),
79 fResidualEta(0x0), fResidualPhi(0x0), fCutEtaPhiSum(kFALSE), fCutEtaPhiSeparate(kFALSE),
80 fCutR(0), fCutEta(0), fCutPhi(0),
81 fClusterWindow(0), fMass(0),
82 fStepSurface(0), fStepCluster(0),
83 fTrackCutsType(0), fCutMinTrackPt(0), fCutMinNClusterTPC(0),
84 fCutMinNClusterITS(0), fCutMaxChi2PerClusterTPC(0), fCutMaxChi2PerClusterITS(0),
85 fCutRequireTPCRefit(kFALSE), fCutRequireITSRefit(kFALSE), fCutAcceptKinkDaughters(kFALSE),
86 fCutMaxDCAToVertexXY(0), fCutMaxDCAToVertexZ(0), fCutDCAToVertex2D(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();
107 //______________________________________________________________________
108 AliEMCALRecoUtils::AliEMCALRecoUtils(const AliEMCALRecoUtils & reco)
110 fParticleType(reco.fParticleType), fPosAlgo(reco.fPosAlgo), fW0(reco.fW0),
111 fNonLinearityFunction(reco.fNonLinearityFunction), fNonLinearThreshold(reco.fNonLinearThreshold),
112 fSmearClusterEnergy(reco.fSmearClusterEnergy), fRandom(),
113 fCellsRecalibrated(reco.fCellsRecalibrated),
114 fRecalibration(reco.fRecalibration), fEMCALRecalibrationFactors(reco.fEMCALRecalibrationFactors),
115 fTimeRecalibration(reco.fTimeRecalibration), fEMCALTimeRecalibrationFactors(reco.fEMCALTimeRecalibrationFactors),
116 fUseRunCorrectionFactors(reco.fUseRunCorrectionFactors), fRunCorrectionFactorsSet(reco.fRunCorrectionFactorsSet),
117 fRemoveBadChannels(reco.fRemoveBadChannels), fRecalDistToBadChannels(reco.fRecalDistToBadChannels),
118 fEMCALBadChannelMap(reco.fEMCALBadChannelMap),
119 fNCellsFromEMCALBorder(reco.fNCellsFromEMCALBorder), fNoEMCALBorderAtEta0(reco.fNoEMCALBorderAtEta0),
120 fRejectExoticCluster(reco.fRejectExoticCluster), fRejectExoticCells(reco.fRejectExoticCells),
121 fExoticCellFraction(reco.fExoticCellFraction), fExoticCellDiffTime(reco.fExoticCellDiffTime),
122 fExoticCellMinAmplitude(reco.fExoticCellMinAmplitude),
123 fPIDUtils(reco.fPIDUtils), fAODFilterMask(reco.fAODFilterMask),
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 fTrackCutsType(reco.fTrackCutsType), fCutMinTrackPt(reco.fCutMinTrackPt),
133 fCutMinNClusterTPC(reco.fCutMinNClusterTPC), fCutMinNClusterITS(reco.fCutMinNClusterITS),
134 fCutMaxChi2PerClusterTPC(reco.fCutMaxChi2PerClusterTPC), fCutMaxChi2PerClusterITS(reco.fCutMaxChi2PerClusterITS),
135 fCutRequireTPCRefit(reco.fCutRequireTPCRefit), fCutRequireITSRefit(reco.fCutRequireITSRefit),
136 fCutAcceptKinkDaughters(reco.fCutAcceptKinkDaughters), fCutMaxDCAToVertexXY(reco.fCutMaxDCAToVertexXY),
137 fCutMaxDCAToVertexZ(reco.fCutMaxDCAToVertexZ), fCutDCAToVertex2D(reco.fCutDCAToVertex2D)
141 for(Int_t i = 0; i < 15 ; i++) { fMisalRotShift[i] = reco.fMisalRotShift[i] ;
142 fMisalTransShift[i] = reco.fMisalTransShift[i] ; }
143 for(Int_t i = 0; i < 7 ; i++) { fNonLinearityParams[i] = reco.fNonLinearityParams[i] ; }
144 for(Int_t i = 0; i < 3 ; i++) { fSmearClusterParam[i] = reco.fSmearClusterParam[i] ; }
149 //______________________________________________________________________
150 AliEMCALRecoUtils & AliEMCALRecoUtils::operator = (const AliEMCALRecoUtils & reco)
152 //Assignment operator
154 if(this == &reco)return *this;
155 ((TNamed *)this)->operator=(reco);
157 for(Int_t i = 0; i < 15 ; i++) { fMisalTransShift[i] = reco.fMisalTransShift[i] ;
158 fMisalRotShift[i] = reco.fMisalRotShift[i] ; }
159 for(Int_t i = 0; i < 7 ; i++) { fNonLinearityParams[i] = reco.fNonLinearityParams[i] ; }
160 for(Int_t i = 0; i < 3 ; i++) { fSmearClusterParam[i] = reco.fSmearClusterParam[i] ; }
162 fParticleType = reco.fParticleType;
163 fPosAlgo = reco.fPosAlgo;
166 fNonLinearityFunction = reco.fNonLinearityFunction;
167 fNonLinearThreshold = reco.fNonLinearThreshold;
168 fSmearClusterEnergy = reco.fSmearClusterEnergy;
170 fCellsRecalibrated = reco.fCellsRecalibrated;
171 fRecalibration = reco.fRecalibration;
172 fEMCALRecalibrationFactors = reco.fEMCALRecalibrationFactors;
174 fTimeRecalibration = reco.fTimeRecalibration;
175 fEMCALTimeRecalibrationFactors = reco.fEMCALTimeRecalibrationFactors;
177 fUseRunCorrectionFactors = reco.fUseRunCorrectionFactors;
178 fRunCorrectionFactorsSet = reco.fRunCorrectionFactorsSet;
180 fRemoveBadChannels = reco.fRemoveBadChannels;
181 fRecalDistToBadChannels = reco.fRecalDistToBadChannels;
182 fEMCALBadChannelMap = reco.fEMCALBadChannelMap;
184 fNCellsFromEMCALBorder = reco.fNCellsFromEMCALBorder;
185 fNoEMCALBorderAtEta0 = reco.fNoEMCALBorderAtEta0;
187 fRejectExoticCluster = reco.fRejectExoticCluster;
188 fRejectExoticCells = reco.fRejectExoticCells;
189 fExoticCellFraction = reco.fExoticCellFraction;
190 fExoticCellDiffTime = reco.fExoticCellDiffTime;
191 fExoticCellMinAmplitude = reco.fExoticCellMinAmplitude;
193 fPIDUtils = reco.fPIDUtils;
195 fAODFilterMask = reco.fAODFilterMask;
197 fCutEtaPhiSum = reco.fCutEtaPhiSum;
198 fCutEtaPhiSeparate = reco.fCutEtaPhiSeparate;
200 fCutEta = reco.fCutEta;
201 fCutPhi = reco.fCutPhi;
202 fClusterWindow = reco.fClusterWindow;
204 fStepSurface = reco.fStepSurface;
205 fStepCluster = reco.fStepCluster;
207 fTrackCutsType = reco.fTrackCutsType;
208 fCutMinTrackPt = reco.fCutMinTrackPt;
209 fCutMinNClusterTPC = reco.fCutMinNClusterTPC;
210 fCutMinNClusterITS = reco.fCutMinNClusterITS;
211 fCutMaxChi2PerClusterTPC = reco.fCutMaxChi2PerClusterTPC;
212 fCutMaxChi2PerClusterITS = reco.fCutMaxChi2PerClusterITS;
213 fCutRequireTPCRefit = reco.fCutRequireTPCRefit;
214 fCutRequireITSRefit = reco.fCutRequireITSRefit;
215 fCutAcceptKinkDaughters = reco.fCutAcceptKinkDaughters;
216 fCutMaxDCAToVertexXY = reco.fCutMaxDCAToVertexXY;
217 fCutMaxDCAToVertexZ = reco.fCutMaxDCAToVertexZ;
218 fCutDCAToVertex2D = reco.fCutDCAToVertex2D;
220 if(reco.fResidualEta)
222 // assign or copy construct
225 *fResidualEta = *reco.fResidualEta;
229 fResidualEta = new TArrayF(*reco.fResidualEta);
234 if(fResidualEta)delete fResidualEta;
238 if(reco.fResidualPhi)
240 // assign or copy construct
243 *fResidualPhi = *reco.fResidualPhi;
247 fResidualPhi = new TArrayF(*reco.fResidualPhi);
252 if(fResidualPhi)delete fResidualPhi;
256 if(reco.fMatchedTrackIndex)
258 // assign or copy construct
259 if(fMatchedTrackIndex)
261 *fMatchedTrackIndex = *reco.fMatchedTrackIndex;
265 fMatchedTrackIndex = new TArrayI(*reco.fMatchedTrackIndex);
270 if(fMatchedTrackIndex)delete fMatchedTrackIndex;
271 fMatchedTrackIndex = 0;
274 if(reco.fMatchedClusterIndex)
276 // assign or copy construct
277 if(fMatchedClusterIndex)
279 *fMatchedClusterIndex = *reco.fMatchedClusterIndex;
283 fMatchedClusterIndex = new TArrayI(*reco.fMatchedClusterIndex);
288 if(fMatchedClusterIndex)delete fMatchedClusterIndex;
289 fMatchedClusterIndex = 0;
296 //_____________________________________
297 AliEMCALRecoUtils::~AliEMCALRecoUtils()
301 if(fEMCALRecalibrationFactors)
303 fEMCALRecalibrationFactors->Clear();
304 delete fEMCALRecalibrationFactors;
307 if(fEMCALTimeRecalibrationFactors)
309 fEMCALTimeRecalibrationFactors->Clear();
310 delete fEMCALTimeRecalibrationFactors;
313 if(fEMCALBadChannelMap)
315 fEMCALBadChannelMap->Clear();
316 delete fEMCALBadChannelMap;
319 delete fMatchedTrackIndex ;
320 delete fMatchedClusterIndex ;
321 delete fResidualEta ;
322 delete fResidualPhi ;
328 //_______________________________________________________________________________
329 Bool_t AliEMCALRecoUtils::AcceptCalibrateCell(const Int_t absID, const Int_t bc,
330 Float_t & amp, Double_t & time,
331 AliVCaloCells* cells)
333 // Reject cell if criteria not passed and calibrate it
335 AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
337 if(absID < 0 || absID >= 24*48*geom->GetNumberOfSuperModules()) return kFALSE;
339 Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1;
341 if(!geom->GetCellIndex(absID,imod,iTower,iIphi,iIeta))
343 // cell absID does not exist
348 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);
350 // Do not include bad channels found in analysis,
351 if( IsBadChannelsRemovalSwitchedOn() && GetEMCALChannelStatus(imod, ieta, iphi))
357 amp = cells->GetCellAmplitude(absID);
358 if(!fCellsRecalibrated && IsRecalibrationOn())
359 amp *= GetEMCALChannelRecalibrationFactor(imod,ieta,iphi);
363 time = cells->GetCellTime(absID);
365 RecalibrateCellTime(absID,bc,time);
370 //_____________________________________________________________________________
371 Bool_t AliEMCALRecoUtils::CheckCellFiducialRegion(const AliEMCALGeometry* geom,
372 const AliVCluster* cluster,
373 AliVCaloCells* cells)
375 // Given the list of AbsId of the cluster, get the maximum cell and
376 // check if there are fNCellsFromBorder from the calorimeter border
380 AliInfo("Cluster pointer null!");
384 //If the distance to the border is 0 or negative just exit accept all clusters
385 if(cells->GetType()==AliVCaloCells::kEMCALCell && fNCellsFromEMCALBorder <= 0 ) return kTRUE;
387 Int_t absIdMax = -1, iSM =-1, ieta = -1, iphi = -1;
388 Bool_t shared = kFALSE;
389 GetMaxEnergyCell(geom, cells, cluster, absIdMax, iSM, ieta, iphi, shared);
391 AliDebug(2,Form("Cluster Max AbsId %d, Cell Energy %2.2f, Cluster Energy %2.2f, Ncells from border %d, EMCAL eta=0 %d\n",
392 absIdMax, cells->GetCellAmplitude(absIdMax), cluster->E(), fNCellsFromEMCALBorder, fNoEMCALBorderAtEta0));
394 if(absIdMax==-1) return kFALSE;
396 //Check if the cell is close to the borders:
397 Bool_t okrow = kFALSE;
398 Bool_t okcol = kFALSE;
400 if(iSM < 0 || iphi < 0 || ieta < 0 )
402 AliFatal(Form("Negative value for super module: %d, or cell ieta: %d, or cell iphi: %d, check EMCAL geometry name\n",
409 if(iphi >= fNCellsFromEMCALBorder && iphi < 24-fNCellsFromEMCALBorder) okrow =kTRUE;
411 else if (iSM >=10 && ( ( geom->GetEMCGeometry()->GetGeoName()).Contains("12SMV1")))
413 if(iphi >= fNCellsFromEMCALBorder && iphi < 8-fNCellsFromEMCALBorder) okrow =kTRUE; //1/3 sm case
417 if(iphi >= fNCellsFromEMCALBorder && iphi < 12-fNCellsFromEMCALBorder) okrow =kTRUE; // half SM case
421 if(!fNoEMCALBorderAtEta0)
423 if(ieta > fNCellsFromEMCALBorder && ieta < 48-fNCellsFromEMCALBorder) okcol =kTRUE;
429 if(ieta >= fNCellsFromEMCALBorder) okcol = kTRUE;
433 if(ieta < 48-fNCellsFromEMCALBorder) okcol = kTRUE;
437 AliDebug(2,Form("EMCAL Cluster in %d cells fiducial volume: ieta %d, iphi %d, SM %d: column? %d, row? %d\nq",
438 fNCellsFromEMCALBorder, ieta, iphi, iSM, okcol, okrow));
442 //printf("Accept\n");
447 //printf("Reject\n");
448 AliDebug(2,Form("Reject cluster in border, max cell : ieta %d, iphi %d, SM %d\n",ieta, iphi, iSM));
455 //_______________________________________________________________________________
456 Bool_t AliEMCALRecoUtils::ClusterContainsBadChannel(const AliEMCALGeometry* geom,
457 const UShort_t* cellList,
460 // Check that in the cluster cells, there is no bad channel of those stored
461 // in fEMCALBadChannelMap or fPHOSBadChannelMap
463 if(!fRemoveBadChannels) return kFALSE;
464 if(!fEMCALBadChannelMap) return kFALSE;
469 for(Int_t iCell = 0; iCell<nCells; iCell++)
471 //Get the column and row
472 Int_t iTower = -1, iIphi = -1, iIeta = -1;
473 geom->GetCellIndex(cellList[iCell],imod,iTower,iIphi,iIeta);
474 if(fEMCALBadChannelMap->GetEntries() <= imod) continue;
475 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,irow,icol);
476 if(GetEMCALChannelStatus(imod, icol, irow))
478 AliDebug(2,Form("Cluster with bad channel: SM %d, col %d, row %d\n",imod, icol, irow));
482 }// cell cluster loop
487 //_____________________________________________________________________________________________
488 Bool_t AliEMCALRecoUtils::IsExoticCell(const Int_t absID, AliVCaloCells* cells, const Int_t bc)
490 // Look to cell neighbourhood and reject if it seems exotic
491 // Do before recalibrating the cells
493 if(!fRejectExoticCells) return kFALSE;
495 AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
497 Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1;
498 geom->GetCellIndex(absID,imod,iTower,iIphi,iIeta);
499 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);
501 //Get close cells index, energy and time, not in corners
503 Int_t absID1 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi+1, ieta);
504 Int_t absID2 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi-1, ieta);
506 // In case of cell in eta = 0 border, depending on SM shift the cross cell index
510 if ( ieta == AliEMCALGeoParams::fgkEMCALCols - 1 && !(imod%2) )
512 absID3 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, 0);
513 absID4 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta-1);
515 else if( ieta == 0 && imod%2 )
517 absID3 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta+1);
518 absID4 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, AliEMCALGeoParams::fgkEMCALCols-1);
522 absID3 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta+1);
523 absID4 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta-1);
526 Float_t ecell = 0, ecell1 = 0, ecell2 = 0, ecell3 = 0, ecell4 = 0;
527 Double_t tcell = 0, tcell1 = 0, tcell2 = 0, tcell3 = 0, tcell4 = 0;
528 Bool_t accept = 0, accept1 = 0, accept2 = 0, accept3 = 0, accept4 = 0;
530 accept = AcceptCalibrateCell(absID, bc, ecell ,tcell ,cells);
532 if(!accept) return kTRUE; // reject this cell
534 if(ecell < fExoticCellMinAmplitude) return kFALSE; // do not reject low energy cells
536 accept1 = AcceptCalibrateCell(absID1,bc, ecell1,tcell1,cells);
537 accept2 = AcceptCalibrateCell(absID2,bc, ecell2,tcell2,cells);
538 accept3 = AcceptCalibrateCell(absID3,bc, ecell3,tcell3,cells);
539 accept4 = AcceptCalibrateCell(absID4,bc, ecell4,tcell4,cells);
542 printf("Cell absID %d \n",absID);
543 printf("\t accept1 %d, accept2 %d, accept3 %d, accept4 %d\n",
544 accept1,accept2,accept3,accept4);
545 printf("\t id %d: id1 %d, id2 %d, id3 %d, id4 %d\n",
546 absID,absID1,absID2,absID3,absID4);
547 printf("\t e %f: e1 %f, e2 %f, e3 %f, e4 %f\n",
548 ecell,ecell1,ecell2,ecell3,ecell4);
549 printf("\t t %f: t1 %f, t2 %f, t3 %f, t4 %f;\n dt1 %f, dt2 %f, dt3 %f, dt4 %f\n",
550 tcell*1.e9,tcell1*1.e9,tcell2*1.e9,tcell3*1.e9,tcell4*1.e9,
551 TMath::Abs(tcell-tcell1)*1.e9, TMath::Abs(tcell-tcell2)*1.e9, TMath::Abs(tcell-tcell3)*1.e9, TMath::Abs(tcell-tcell4)*1.e9);
554 if(TMath::Abs(tcell-tcell1)*1.e9 > fExoticCellDiffTime) ecell1 = 0 ;
555 if(TMath::Abs(tcell-tcell2)*1.e9 > fExoticCellDiffTime) ecell2 = 0 ;
556 if(TMath::Abs(tcell-tcell3)*1.e9 > fExoticCellDiffTime) ecell3 = 0 ;
557 if(TMath::Abs(tcell-tcell4)*1.e9 > fExoticCellDiffTime) ecell4 = 0 ;
559 Float_t eCross = ecell1+ecell2+ecell3+ecell4;
561 //printf("\t eCell %f, eCross %f, 1-eCross/eCell %f\n",ecell,eCross,1-eCross/ecell);
563 if(1-eCross/ecell > fExoticCellFraction)
565 AliDebug(2,Form("AliEMCALRecoUtils::IsExoticCell() - EXOTIC CELL id %d, eCell %f, eCross %f, 1-eCross/eCell %f\n",
566 absID,ecell,eCross,1-eCross/ecell));
573 //___________________________________________________________________
574 Bool_t AliEMCALRecoUtils::IsExoticCluster(const AliVCluster *cluster,
575 AliVCaloCells *cells,
578 // Check if the cluster highest energy tower is exotic
582 AliInfo("Cluster pointer null!");
586 if(!fRejectExoticCluster) return kFALSE;
588 // Get highest energy tower
589 AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
590 Int_t iSupMod = -1, absId = -1, ieta = -1, iphi = -1;
591 Bool_t shared = kFALSE;
592 GetMaxEnergyCell(geom, cells, cluster, absId, iSupMod, ieta, iphi, shared);
594 return IsExoticCell(absId,cells,bc);
598 //_______________________________________________________________________
599 Float_t AliEMCALRecoUtils::SmearClusterEnergy(const AliVCluster* cluster)
601 //In case of MC analysis, smear energy to match resolution/calibration in real data
605 AliInfo("Cluster pointer null!");
609 Float_t energy = cluster->E() ;
610 Float_t rdmEnergy = energy ;
611 if(fSmearClusterEnergy)
613 rdmEnergy = fRandom.Gaus(energy,fSmearClusterParam[0] * TMath::Sqrt(energy) +
614 fSmearClusterParam[1] * energy +
615 fSmearClusterParam[2] );
616 AliDebug(2, Form("Energy: original %f, smeared %f\n", energy, rdmEnergy));
622 //____________________________________________________________________________
623 Float_t AliEMCALRecoUtils::CorrectClusterEnergyLinearity(AliVCluster* cluster)
625 // Correct cluster energy from non linearity functions
629 AliInfo("Cluster pointer null!");
633 Float_t energy = cluster->E();
635 switch (fNonLinearityFunction)
640 //Non-Linearity correction (from MC with function ([0]*exp(-[1]/E))+(([2]/([3]*2.*TMath::Pi())*exp(-(E-[4])^2/(2.*[3]^2)))))
641 //Double_t fNonLinearityParams[0] = 1.014;
642 //Double_t fNonLinearityParams[1] = -0.03329;
643 //Double_t fNonLinearityParams[2] = -0.3853;
644 //Double_t fNonLinearityParams[3] = 0.5423;
645 //Double_t fNonLinearityParams[4] = -0.4335;
646 energy *= (fNonLinearityParams[0]*exp(-fNonLinearityParams[1]/energy))+
647 ((fNonLinearityParams[2]/(fNonLinearityParams[3]*2.*TMath::Pi())*
648 exp(-(energy-fNonLinearityParams[4])*(energy-fNonLinearityParams[4])/(2.*fNonLinearityParams[3]*fNonLinearityParams[3]))));
654 //Non-Linearity correction (from Olga Data with function p0+p1*exp(-p2*E))
655 //Double_t fNonLinearityParams[0] = 1.04;
656 //Double_t fNonLinearityParams[1] = -0.1445;
657 //Double_t fNonLinearityParams[2] = 1.046;
658 energy /= (fNonLinearityParams[0]+fNonLinearityParams[1]*exp(-fNonLinearityParams[2]*energy)); //Olga function
662 case kPi0GammaConversion:
664 //Non-Linearity correction (Nicolas from Dimitri Data with function C*[1-a*exp(-b*E)])
665 //fNonLinearityParams[0] = 0.139393/0.1349766;
666 //fNonLinearityParams[1] = 0.0566186;
667 //fNonLinearityParams[2] = 0.982133;
668 energy /= fNonLinearityParams[0]*(1-fNonLinearityParams[1]*exp(-fNonLinearityParams[2]*energy));
675 //From beam test, Alexei's results, for different ZS thresholds
676 // th=30 MeV; th = 45 MeV; th = 75 MeV
677 //fNonLinearityParams[0] = 1.007; 1.003; 1.002
678 //fNonLinearityParams[1] = 0.894; 0.719; 0.797
679 //fNonLinearityParams[2] = 0.246; 0.334; 0.358
680 //Rescale the param[0] with 1.03
681 energy /= fNonLinearityParams[0]/(1+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]));
686 case kBeamTestCorrected:
688 //From beam test, corrected for material between beam and EMCAL
689 //fNonLinearityParams[0] = 0.99078
690 //fNonLinearityParams[1] = 0.161499;
691 //fNonLinearityParams[2] = 0.655166;
692 //fNonLinearityParams[3] = 0.134101;
693 //fNonLinearityParams[4] = 163.282;
694 //fNonLinearityParams[5] = 23.6904;
695 //fNonLinearityParams[6] = 0.978;
696 energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
702 AliDebug(2,"No correction on the energy\n");
710 //__________________________________________________
711 void AliEMCALRecoUtils::InitNonLinearityParam()
713 //Initialising Non Linearity Parameters
715 if(fNonLinearityFunction == kPi0MC)
717 fNonLinearityParams[0] = 1.014;
718 fNonLinearityParams[1] = -0.03329;
719 fNonLinearityParams[2] = -0.3853;
720 fNonLinearityParams[3] = 0.5423;
721 fNonLinearityParams[4] = -0.4335;
724 if(fNonLinearityFunction == kPi0GammaGamma)
726 fNonLinearityParams[0] = 1.04;
727 fNonLinearityParams[1] = -0.1445;
728 fNonLinearityParams[2] = 1.046;
731 if(fNonLinearityFunction == kPi0GammaConversion)
733 fNonLinearityParams[0] = 0.139393;
734 fNonLinearityParams[1] = 0.0566186;
735 fNonLinearityParams[2] = 0.982133;
738 if(fNonLinearityFunction == kBeamTest)
740 if(fNonLinearThreshold == 30)
742 fNonLinearityParams[0] = 1.007;
743 fNonLinearityParams[1] = 0.894;
744 fNonLinearityParams[2] = 0.246;
746 if(fNonLinearThreshold == 45)
748 fNonLinearityParams[0] = 1.003;
749 fNonLinearityParams[1] = 0.719;
750 fNonLinearityParams[2] = 0.334;
752 if(fNonLinearThreshold == 75)
754 fNonLinearityParams[0] = 1.002;
755 fNonLinearityParams[1] = 0.797;
756 fNonLinearityParams[2] = 0.358;
760 if(fNonLinearityFunction == kBeamTestCorrected)
762 fNonLinearityParams[0] = 0.99078;
763 fNonLinearityParams[1] = 0.161499;
764 fNonLinearityParams[2] = 0.655166;
765 fNonLinearityParams[3] = 0.134101;
766 fNonLinearityParams[4] = 163.282;
767 fNonLinearityParams[5] = 23.6904;
768 fNonLinearityParams[6] = 0.978;
772 //_________________________________________________________
773 Float_t AliEMCALRecoUtils::GetDepth(const Float_t energy,
774 const Int_t iParticle,
775 const Int_t iSM) const
777 //Calculate shower depth for a given cluster energy and particle type
787 depth = x0 * (TMath::Log(energy*1000/ ecr) + 0.5); //Multiply energy by 1000 to transform to MeV
791 depth = x0 * (TMath::Log(energy*1000/ ecr) - 0.5); //Multiply energy by 1000 to transform to MeV
799 gGeoManager->cd("ALIC_1/XEN1_1");
800 TGeoNode *geoXEn1 = gGeoManager->GetCurrentNode();
801 TGeoNodeMatrix *geoSM = dynamic_cast<TGeoNodeMatrix *>(geoXEn1->GetDaughter(iSM));
804 TGeoVolume *geoSMVol = geoSM->GetVolume();
805 TGeoShape *geoSMShape = geoSMVol->GetShape();
806 TGeoBBox *geoBox = dynamic_cast<TGeoBBox *>(geoSMShape);
807 if(geoBox) depth = 0.5 * geoBox->GetDX()*2 ;
808 else AliFatal("Null GEANT box");
810 else AliFatal("NULL GEANT node matrix");
814 depth = x0 * (TMath::Log(energy*1000 / ecr) - 0.5); //Multiply energy by 1000 to transform to MeV
820 depth = x0 * (TMath::Log(energy*1000 / ecr) + 0.5); //Multiply energy by 1000 to transform to MeV
826 //____________________________________________________________________
827 void AliEMCALRecoUtils::GetMaxEnergyCell(const AliEMCALGeometry *geom,
828 AliVCaloCells* cells,
829 const AliVCluster* clu,
836 //For a given CaloCluster gets the absId of the cell
837 //with maximum energy deposit.
840 Double_t eCell = -1.;
841 Float_t fraction = 1.;
842 Float_t recalFactor = 1.;
843 Int_t cellAbsId = -1 ;
852 AliInfo("Cluster pointer null!");
853 absId=-1; iSupMod0=-1, ieta = -1; iphi = -1; shared = -1;
857 for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++)
859 cellAbsId = clu->GetCellAbsId(iDig);
860 fraction = clu->GetCellAmplitudeFraction(iDig);
861 //printf("a Cell %d, id, %d, amp %f, fraction %f\n",iDig,cellAbsId,cells->GetCellAmplitude(cellAbsId),fraction);
862 if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
863 geom->GetCellIndex(cellAbsId,iSupMod,iTower,iIphi,iIeta);
864 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);
869 else if(iSupMod0!=iSupMod)
872 //printf("AliEMCALRecoUtils::GetMaxEnergyCell() - SHARED CLUSTER\n");
874 if(!fCellsRecalibrated && IsRecalibrationOn())
876 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
878 eCell = cells->GetCellAmplitude(cellAbsId)*fraction*recalFactor;
879 //printf("b Cell %d, id, %d, amp %f, fraction %f\n",iDig,cellAbsId,eCell,fraction);
884 //printf("\t new max: cell %d, e %f, ecell %f\n",maxId, eMax,eCell);
888 //Get from the absid the supermodule, tower and eta/phi numbers
889 geom->GetCellIndex(absId,iSupMod,iTower,iIphi,iIeta);
890 //Gives SuperModule and Tower numbers
891 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,
892 iIphi, iIeta,iphi,ieta);
893 //printf("Max id %d, iSM %d, col %d, row %d\n",absId,iSupMod,ieta,iphi);
894 //printf("Max end---\n");
897 //______________________________________
898 void AliEMCALRecoUtils::InitParameters()
900 // Initialize data members with default values
902 fParticleType = kPhoton;
903 fPosAlgo = kUnchanged;
906 fNonLinearityFunction = kNoCorrection;
907 fNonLinearThreshold = 30;
909 fExoticCellFraction = 0.97;
910 fExoticCellDiffTime = 1e6;
911 fExoticCellMinAmplitude = 0.5;
915 fCutEtaPhiSum = kTRUE;
916 fCutEtaPhiSeparate = kFALSE;
922 fClusterWindow = 100;
927 fTrackCutsType = kLooseCut;
930 fCutMinNClusterTPC = -1;
931 fCutMinNClusterITS = -1;
933 fCutMaxChi2PerClusterTPC = 1e10;
934 fCutMaxChi2PerClusterITS = 1e10;
936 fCutRequireTPCRefit = kFALSE;
937 fCutRequireITSRefit = kFALSE;
938 fCutAcceptKinkDaughters = kFALSE;
940 fCutMaxDCAToVertexXY = 1e10;
941 fCutMaxDCAToVertexZ = 1e10;
942 fCutDCAToVertex2D = kFALSE;
945 //Misalignment matrices
946 for(Int_t i = 0; i < 15 ; i++)
948 fMisalTransShift[i] = 0.;
949 fMisalRotShift[i] = 0.;
953 for(Int_t i = 0; i < 7 ; i++) fNonLinearityParams[i] = 0.;
955 //For kBeamTestCorrected case, but default is no correction
956 fNonLinearityParams[0] = 0.99078;
957 fNonLinearityParams[1] = 0.161499;
958 fNonLinearityParams[2] = 0.655166;
959 fNonLinearityParams[3] = 0.134101;
960 fNonLinearityParams[4] = 163.282;
961 fNonLinearityParams[5] = 23.6904;
962 fNonLinearityParams[6] = 0.978;
964 //For kPi0GammaGamma case
965 //fNonLinearityParams[0] = 0.1457/0.1349766/1.038;
966 //fNonLinearityParams[1] = -0.02024/0.1349766/1.038;
967 //fNonLinearityParams[2] = 1.046;
969 //Cluster energy smearing
970 fSmearClusterEnergy = kFALSE;
971 fSmearClusterParam[0] = 0.07; // * sqrt E term
972 fSmearClusterParam[1] = 0.00; // * E term
973 fSmearClusterParam[2] = 0.00; // constant
976 //_____________________________________________________
977 void AliEMCALRecoUtils::InitEMCALRecalibrationFactors()
979 //Init EMCAL recalibration factors
980 AliDebug(2,"AliCalorimeterUtils::InitEMCALRecalibrationFactors()");
981 //In order to avoid rewriting the same histograms
982 Bool_t oldStatus = TH1::AddDirectoryStatus();
983 TH1::AddDirectory(kFALSE);
985 fEMCALRecalibrationFactors = new TObjArray(12);
986 for (int i = 0; i < 12; i++)
987 fEMCALRecalibrationFactors->Add(new TH2F(Form("EMCALRecalFactors_SM%d",i),
988 Form("EMCALRecalFactors_SM%d",i), 48, 0, 48, 24, 0, 24));
989 //Init the histograms with 1
990 for (Int_t sm = 0; sm < 12; sm++)
992 for (Int_t i = 0; i < 48; i++)
994 for (Int_t j = 0; j < 24; j++)
996 SetEMCALChannelRecalibrationFactor(sm,i,j,1.);
1001 fEMCALRecalibrationFactors->SetOwner(kTRUE);
1002 fEMCALRecalibrationFactors->Compress();
1004 //In order to avoid rewriting the same histograms
1005 TH1::AddDirectory(oldStatus);
1008 //_________________________________________________________
1009 void AliEMCALRecoUtils::InitEMCALTimeRecalibrationFactors()
1011 //Init EMCAL recalibration factors
1012 AliDebug(2,"AliCalorimeterUtils::InitEMCALRecalibrationFactors()");
1013 //In order to avoid rewriting the same histograms
1014 Bool_t oldStatus = TH1::AddDirectoryStatus();
1015 TH1::AddDirectory(kFALSE);
1017 fEMCALTimeRecalibrationFactors = new TObjArray(4);
1018 for (int i = 0; i < 4; i++)
1019 fEMCALTimeRecalibrationFactors->Add(new TH1F(Form("hAllTimeAvBC%d",i),
1020 Form("hAllTimeAvBC%d",i),
1021 48*24*12,0.,48*24*12) );
1022 //Init the histograms with 1
1023 for (Int_t bc = 0; bc < 4; bc++)
1025 for (Int_t i = 0; i < 48*24*12; i++)
1026 SetEMCALChannelTimeRecalibrationFactor(bc,i,0.);
1029 fEMCALTimeRecalibrationFactors->SetOwner(kTRUE);
1030 fEMCALTimeRecalibrationFactors->Compress();
1032 //In order to avoid rewriting the same histograms
1033 TH1::AddDirectory(oldStatus);
1036 //____________________________________________________
1037 void AliEMCALRecoUtils::InitEMCALBadChannelStatusMap()
1039 //Init EMCAL bad channels map
1040 AliDebug(2,"AliEMCALRecoUtils::InitEMCALBadChannelStatusMap()");
1041 //In order to avoid rewriting the same histograms
1042 Bool_t oldStatus = TH1::AddDirectoryStatus();
1043 TH1::AddDirectory(kFALSE);
1045 fEMCALBadChannelMap = new TObjArray(12);
1046 //TH2F * hTemp = new TH2I("EMCALBadChannelMap","EMCAL SuperModule bad channel map", 48, 0, 48, 24, 0, 24);
1047 for (int i = 0; i < 12; i++)
1049 fEMCALBadChannelMap->Add(new TH2I(Form("EMCALBadChannelMap_Mod%d",i),Form("EMCALBadChannelMap_Mod%d",i), 48, 0, 48, 24, 0, 24));
1052 fEMCALBadChannelMap->SetOwner(kTRUE);
1053 fEMCALBadChannelMap->Compress();
1055 //In order to avoid rewriting the same histograms
1056 TH1::AddDirectory(oldStatus);
1059 //____________________________________________________________________________
1060 void AliEMCALRecoUtils::RecalibrateClusterEnergy(const AliEMCALGeometry* geom,
1061 AliVCluster * cluster,
1062 AliVCaloCells * cells,
1065 // Recalibrate the cluster energy and Time, considering the recalibration map
1066 // and the energy of the cells and time that compose the cluster.
1067 // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1071 AliInfo("Cluster pointer null!");
1075 //Get the cluster number of cells and list of absId, check what kind of cluster do we have.
1076 UShort_t * index = cluster->GetCellsAbsId() ;
1077 Double_t * fraction = cluster->GetCellsAmplitudeFraction() ;
1078 Int_t ncells = cluster->GetNCells();
1080 //Initialize some used variables
1083 Int_t icol =-1, irow =-1, imod=1;
1084 Float_t factor = 1, frac = 0;
1085 Int_t absIdMax = -1;
1088 //Loop on the cells, get the cell amplitude and recalibration factor, multiply and and to the new energy
1089 for(Int_t icell = 0; icell < ncells; icell++)
1091 absId = index[icell];
1092 frac = fraction[icell];
1093 if(frac < 1e-5) frac = 1; //in case of EMCAL, this is set as 0 since unfolding is off
1095 if(!fCellsRecalibrated && IsRecalibrationOn())
1098 Int_t iTower = -1, iIphi = -1, iIeta = -1;
1099 geom->GetCellIndex(absId,imod,iTower,iIphi,iIeta);
1100 if(fEMCALRecalibrationFactors->GetEntries() <= imod) continue;
1101 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,irow,icol);
1102 factor = GetEMCALChannelRecalibrationFactor(imod,icol,irow);
1104 AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - recalibrate cell: module %d, col %d, row %d, cell fraction %f,recalibration factor %f, cell energy %f\n",
1105 imod,icol,irow,frac,factor,cells->GetCellAmplitude(absId)));
1109 energy += cells->GetCellAmplitude(absId)*factor*frac;
1111 if(emax < cells->GetCellAmplitude(absId)*factor*frac)
1113 emax = cells->GetCellAmplitude(absId)*factor*frac;
1118 cluster->SetE(energy);
1120 AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - Energy before %f, after %f\n",cluster->E(),energy));
1122 // Recalculate time of cluster only for ESDs
1123 if(!strcmp("AliESDCaloCluster",Form("%s",cluster->ClassName())))
1126 Double_t weightedTime = 0;
1127 Double_t weight = 0;
1128 Double_t weightTot = 0;
1129 Double_t maxcellTime = 0;
1130 for(Int_t icell = 0; icell < ncells; icell++)
1132 absId = index[icell];
1133 frac = fraction[icell];
1134 if(frac < 1e-5) frac = 1; //in case of EMCAL, this is set as 0 since unfolding is off
1136 Double_t celltime = cells->GetCellTime(absId);
1137 RecalibrateCellTime(absId, bc, celltime);
1138 if(absId == absIdMax) maxcellTime = celltime;
1140 if(!fCellsRecalibrated)
1142 Int_t iTower = -1, iIphi = -1, iIeta = -1;
1143 geom->GetCellIndex(absId,imod,iTower,iIphi,iIeta);
1144 if(fEMCALRecalibrationFactors->GetEntries() <= imod) continue;
1145 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,irow,icol);
1146 factor = GetEMCALChannelRecalibrationFactor(imod,icol,irow);
1148 AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - recalibrate cell:"
1149 " module %d, col %d, row %d, cell fraction %f,recalibration factor %f, cell energy %f\n",
1150 imod,icol,irow,frac,factor,cells->GetCellTime(absId)));
1154 weight = GetCellWeight(cells->GetCellAmplitude(absId)*factor*frac , energy );
1155 weightTot += weight;
1156 weightedTime += celltime * weight;
1160 cluster->SetTOF(weightedTime/weightTot);
1162 cluster->SetTOF(maxcellTime);
1166 //_____________________________________________________________
1167 void AliEMCALRecoUtils::RecalibrateCells(AliVCaloCells * cells,
1170 // Recalibrate the cells time and energy, considering the recalibration map and the energy
1171 // of the cells that compose the cluster.
1172 // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1174 if(!IsRecalibrationOn() && !IsTimeRecalibrationOn()) return;
1178 AliInfo("Cells pointer null!");
1183 Bool_t accept = kFALSE;
1187 Int_t nEMcell = cells->GetNumberOfCells() ;
1188 for (Int_t iCell = 0; iCell < nEMcell; iCell++)
1190 absId = cells->GetCellNumber(iCell);
1192 accept = AcceptCalibrateCell(absId, bc, ecell ,tcell ,cells);
1200 cells->SetCell(iCell,absId,ecell, tcell);
1203 fCellsRecalibrated = kTRUE;
1206 //_______________________________________________________________________________________________________
1207 void AliEMCALRecoUtils::RecalibrateCellTime(const Int_t absId, const Int_t bc, Double_t & celltime) const
1209 // Recalibrate time of cell with absID considering the recalibration map
1210 // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1212 if(!fCellsRecalibrated && IsTimeRecalibrationOn() && bc >= 0)
1214 celltime -= GetEMCALChannelTimeRecalibrationFactor(bc%4,absId)*1.e-9; ;
1218 //______________________________________________________________________________
1219 void AliEMCALRecoUtils::RecalculateClusterPosition(const AliEMCALGeometry *geom,
1220 AliVCaloCells* cells,
1223 //For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1227 AliInfo("Cluster pointer null!");
1231 if (fPosAlgo==kPosTowerGlobal) RecalculateClusterPositionFromTowerGlobal( geom, cells, clu);
1232 else if(fPosAlgo==kPosTowerIndex) RecalculateClusterPositionFromTowerIndex ( geom, cells, clu);
1233 else AliDebug(2,"Algorithm to recalculate position not selected, do nothing.");
1236 //_____________________________________________________________________________________________
1237 void AliEMCALRecoUtils::RecalculateClusterPositionFromTowerGlobal(const AliEMCALGeometry *geom,
1238 AliVCaloCells* cells,
1241 // For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1242 // The algorithm is a copy of what is done in AliEMCALRecPoint
1244 Double_t eCell = 0.;
1245 Float_t fraction = 1.;
1246 Float_t recalFactor = 1.;
1249 Int_t iTower = -1, iIphi = -1, iIeta = -1;
1250 Int_t iSupModMax = -1, iSM=-1, iphi = -1, ieta = -1;
1251 Float_t weight = 0., totalWeight=0.;
1252 Float_t newPos[3] = {0,0,0};
1253 Double_t pLocal[3], pGlobal[3];
1254 Bool_t shared = kFALSE;
1256 Float_t clEnergy = clu->E(); //Energy already recalibrated previously
1257 GetMaxEnergyCell(geom, cells, clu, absId, iSupModMax, ieta, iphi,shared);
1258 Double_t depth = GetDepth(clEnergy,fParticleType,iSupModMax) ;
1260 //printf("** Cluster energy %f, ncells %d, depth %f\n",clEnergy,clu->GetNCells(),depth);
1262 for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++)
1264 absId = clu->GetCellAbsId(iDig);
1265 fraction = clu->GetCellAmplitudeFraction(iDig);
1266 if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1268 if (!fCellsRecalibrated)
1270 geom->GetCellIndex(absId,iSM,iTower,iIphi,iIeta);
1271 geom->GetCellPhiEtaIndexInSModule(iSM,iTower,iIphi, iIeta,iphi,ieta);
1273 if(IsRecalibrationOn())
1275 recalFactor = GetEMCALChannelRecalibrationFactor(iSM,ieta,iphi);
1279 eCell = cells->GetCellAmplitude(absId)*fraction*recalFactor;
1281 weight = GetCellWeight(eCell,clEnergy);
1282 totalWeight += weight;
1284 geom->RelPosCellInSModule(absId,depth,pLocal[0],pLocal[1],pLocal[2]);
1285 //printf("pLocal (%f,%f,%f), SM %d, absId %d\n",pLocal[0],pLocal[1],pLocal[2],iSupModMax,absId);
1286 geom->GetGlobal(pLocal,pGlobal,iSupModMax);
1287 //printf("pLocal (%f,%f,%f)\n",pGlobal[0],pGlobal[1],pGlobal[2]);
1289 for(int i=0; i<3; i++ ) newPos[i] += (weight*pGlobal[i]);
1294 for(int i=0; i<3; i++ ) newPos[i] /= totalWeight;
1297 //Float_t pos[]={0,0,0};
1298 //clu->GetPosition(pos);
1299 //printf("OldPos : %2.3f,%2.3f,%2.3f\n",pos[0],pos[1],pos[2]);
1300 //printf("NewPos : %2.3f,%2.3f,%2.3f\n",newPos[0],newPos[1],newPos[2]);
1302 if(iSupModMax > 1) //sector 1
1304 newPos[0] +=fMisalTransShift[3];//-=3.093;
1305 newPos[1] +=fMisalTransShift[4];//+=6.82;
1306 newPos[2] +=fMisalTransShift[5];//+=1.635;
1307 //printf(" + : %2.3f,%2.3f,%2.3f\n",fMisalTransShift[3],fMisalTransShift[4],fMisalTransShift[5]);
1310 newPos[0] +=fMisalTransShift[0];//+=1.134;
1311 newPos[1] +=fMisalTransShift[1];//+=8.2;
1312 newPos[2] +=fMisalTransShift[2];//+=1.197;
1313 //printf(" + : %2.3f,%2.3f,%2.3f\n",fMisalTransShift[0],fMisalTransShift[1],fMisalTransShift[2]);
1315 //printf("NewPos : %2.3f,%2.3f,%2.3f\n",newPos[0],newPos[1],newPos[2]);
1317 clu->SetPosition(newPos);
1320 //____________________________________________________________________________________________
1321 void AliEMCALRecoUtils::RecalculateClusterPositionFromTowerIndex(const AliEMCALGeometry *geom,
1322 AliVCaloCells* cells,
1325 // For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1326 // The algorithm works with the tower indeces, averages the indeces and from them it calculates the global position
1328 Double_t eCell = 1.;
1329 Float_t fraction = 1.;
1330 Float_t recalFactor = 1.;
1334 Int_t iIphi = -1, iIeta = -1;
1335 Int_t iSupMod = -1, iSupModMax = -1;
1336 Int_t iphi = -1, ieta =-1;
1337 Bool_t shared = kFALSE;
1339 Float_t clEnergy = clu->E(); //Energy already recalibrated previously.
1340 GetMaxEnergyCell(geom, cells, clu, absId, iSupModMax, ieta, iphi,shared);
1341 Float_t depth = GetDepth(clEnergy,fParticleType,iSupMod) ;
1343 Float_t weight = 0., weightedCol = 0., weightedRow = 0., totalWeight=0.;
1344 Bool_t areInSameSM = kTRUE; //exclude clusters with cells in different SMs for now
1345 Int_t startingSM = -1;
1347 for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++)
1349 absId = clu->GetCellAbsId(iDig);
1350 fraction = clu->GetCellAmplitudeFraction(iDig);
1351 if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1353 if (iDig==0) startingSM = iSupMod;
1354 else if(iSupMod != startingSM) areInSameSM = kFALSE;
1356 eCell = cells->GetCellAmplitude(absId);
1358 geom->GetCellIndex(absId,iSupMod,iTower,iIphi,iIeta);
1359 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);
1361 if (!fCellsRecalibrated)
1363 if(IsRecalibrationOn())
1365 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1369 eCell = cells->GetCellAmplitude(absId)*fraction*recalFactor;
1371 weight = GetCellWeight(eCell,clEnergy);
1372 if(weight < 0) weight = 0;
1373 totalWeight += weight;
1374 weightedCol += ieta*weight;
1375 weightedRow += iphi*weight;
1377 //printf("Max cell? cell %d, amplitude org %f, fraction %f, recalibration %f, amplitude new %f \n",cellAbsId, cells->GetCellAmplitude(cellAbsId), fraction, recalFactor, eCell) ;
1380 Float_t xyzNew[]={0.,0.,0.};
1381 if(areInSameSM == kTRUE)
1383 //printf("In Same SM\n");
1384 weightedCol = weightedCol/totalWeight;
1385 weightedRow = weightedRow/totalWeight;
1386 geom->RecalculateTowerPosition(weightedRow, weightedCol, iSupModMax, depth, fMisalTransShift, fMisalRotShift, xyzNew);
1390 //printf("In Different SM\n");
1391 geom->RecalculateTowerPosition(iphi, ieta, iSupModMax, depth, fMisalTransShift, fMisalRotShift, xyzNew);
1394 clu->SetPosition(xyzNew);
1397 //___________________________________________________________________________________________
1398 void AliEMCALRecoUtils::RecalculateClusterDistanceToBadChannel(const AliEMCALGeometry * geom,
1399 AliVCaloCells* cells,
1400 AliVCluster * cluster)
1402 //re-evaluate distance to bad channel with updated bad map
1404 if(!fRecalDistToBadChannels) return;
1408 AliInfo("Cluster pointer null!");
1412 //Get channels map of the supermodule where the cluster is.
1413 Int_t absIdMax = -1, iSupMod =-1, icolM = -1, irowM = -1;
1414 Bool_t shared = kFALSE;
1415 GetMaxEnergyCell(geom, cells, cluster, absIdMax, iSupMod, icolM, irowM, shared);
1416 TH2D* hMap = (TH2D*)fEMCALBadChannelMap->At(iSupMod);
1419 Float_t minDist = 10000.;
1422 //Loop on tower status map
1423 for(Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++)
1425 for(Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++)
1427 //Check if tower is bad.
1428 if(hMap->GetBinContent(icol,irow)==0) continue;
1429 //printf("AliEMCALRecoUtils::RecalculateDistanceToBadChannels() - \n \t Bad channel in SM %d, col %d, row %d, \n \t Cluster max in col %d, row %d\n",
1430 // iSupMod,icol, irow, icolM,irowM);
1432 dRrow=TMath::Abs(irowM-irow);
1433 dRcol=TMath::Abs(icolM-icol);
1434 dist=TMath::Sqrt(dRrow*dRrow+dRcol*dRcol);
1437 //printf("MIN DISTANCE TO BAD %2.2f\n",dist);
1443 //In case the cluster is shared by 2 SuperModules, need to check the map of the second Super Module
1447 Int_t iSupMod2 = -1;
1449 //The only possible combinations are (0,1), (2,3) ... (8,9)
1450 if(iSupMod%2) iSupMod2 = iSupMod-1;
1451 else iSupMod2 = iSupMod+1;
1452 hMap2 = (TH2D*)fEMCALBadChannelMap->At(iSupMod2);
1454 //Loop on tower status map of second super module
1455 for(Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++)
1457 for(Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++)
1459 //Check if tower is bad.
1460 if(hMap2->GetBinContent(icol,irow)==0) continue;
1461 //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",
1462 // iSupMod2,icol, irow,iSupMod,icolM,irowM);
1463 dRrow=TMath::Abs(irow-irowM);
1467 dRcol=TMath::Abs(icol-(AliEMCALGeoParams::fgkEMCALCols+icolM));
1470 dRcol=TMath::Abs(AliEMCALGeoParams::fgkEMCALCols+icol-icolM);
1473 dist=TMath::Sqrt(dRrow*dRrow+dRcol*dRcol);
1474 if(dist < minDist) minDist = dist;
1477 }// shared cluster in 2 SuperModules
1479 AliDebug(2,Form("Max cluster cell (SM,col,row)=(%d %d %d) - Distance to Bad Channel %2.2f",iSupMod, icolM, irowM, minDist));
1480 cluster->SetDistanceToBadChannel(minDist);
1483 //__________________________________________________________________
1484 void AliEMCALRecoUtils::RecalculateClusterPID(AliVCluster * cluster)
1486 //re-evaluate identification parameters with bayesian
1490 AliInfo("Cluster pointer null!");
1494 if ( cluster->GetM02() != 0)
1495 fPIDUtils->ComputePID(cluster->E(),cluster->GetM02());
1497 Float_t pidlist[AliPID::kSPECIESN+1];
1498 for(Int_t i = 0; i < AliPID::kSPECIESN+1; i++) pidlist[i] = fPIDUtils->GetPIDFinal(i);
1500 cluster->SetPID(pidlist);
1503 //____________________________________________________________________________________________
1504 void AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(const AliEMCALGeometry * geom,
1505 AliVCaloCells* cells,
1506 AliVCluster * cluster)
1508 // Calculates new center of gravity in the local EMCAL-module coordinates
1509 // and tranfers into global ALICE coordinates
1510 // Calculates Dispersion and main axis
1514 AliInfo("Cluster pointer null!");
1520 Double_t eCell = 0.;
1521 Float_t fraction = 1.;
1522 Float_t recalFactor = 1.;
1530 Double_t etai = -1.;
1531 Double_t phii = -1.;
1538 Double_t xmean = 0.;
1539 Double_t zmean = 0.;
1542 for(Int_t iDigit=0; iDigit < cluster->GetNCells(); iDigit++)
1544 //Get from the absid the supermodule, tower and eta/phi numbers
1545 geom->GetCellIndex(cluster->GetCellAbsId(iDigit),iSupMod,iTower,iIphi,iIeta);
1546 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi,iIeta, iphi,ieta);
1548 //Get the cell energy, if recalibration is on, apply factors
1549 fraction = cluster->GetCellAmplitudeFraction(iDigit);
1550 if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1552 if (!fCellsRecalibrated)
1554 if(IsRecalibrationOn())
1556 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1560 eCell = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1562 if(cluster->E() > 0 && eCell > 0)
1564 w = GetCellWeight(eCell,cluster->E());
1566 etai=(Double_t)ieta;
1567 phii=(Double_t)iphi;
1573 dxx += w * etai * etai ;
1575 dzz += w * phii * phii ;
1577 dxz += w * etai * phii ;
1581 AliError(Form("Wrong energy %f and/or amplitude %f\n", eCell, cluster->E()));
1584 //Normalize to the weight
1591 AliError(Form("Wrong weight %f\n", wtot));
1593 //Calculate dispersion
1594 for(Int_t iDigit=0; iDigit < cluster->GetNCells(); iDigit++)
1596 //Get from the absid the supermodule, tower and eta/phi numbers
1597 geom->GetCellIndex(cluster->GetCellAbsId(iDigit),iSupMod,iTower,iIphi,iIeta);
1598 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi,iIeta, iphi,ieta);
1600 //Get the cell energy, if recalibration is on, apply factors
1601 fraction = cluster->GetCellAmplitudeFraction(iDigit);
1602 if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1603 if (IsRecalibrationOn())
1605 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1607 eCell = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1609 if(cluster->E() > 0 && eCell > 0)
1611 w = GetCellWeight(eCell,cluster->E());
1613 etai=(Double_t)ieta;
1614 phii=(Double_t)iphi;
1615 if(w > 0.0) d += w*((etai-xmean)*(etai-xmean)+(phii-zmean)*(phii-zmean));
1618 AliError(Form("Wrong energy %f and/or amplitude %f\n", eCell, cluster->E()));
1621 //Normalize to the weigth and set shower shape parameters
1622 if (wtot > 0 && nstat > 1)
1628 dxx -= xmean * xmean ;
1629 dzz -= zmean * zmean ;
1630 dxz -= xmean * zmean ;
1631 cluster->SetM02(0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ));
1632 cluster->SetM20(0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ));
1637 cluster->SetM20(0.) ;
1638 cluster->SetM02(0.) ;
1642 cluster->SetDispersion(TMath::Sqrt(d)) ;
1644 cluster->SetDispersion(0) ;
1647 //____________________________________________________________________________
1648 void AliEMCALRecoUtils::FindMatches(AliVEvent *event,
1649 TObjArray * clusterArr,
1650 const AliEMCALGeometry *geom)
1652 //This function should be called before the cluster loop
1653 //Before call this function, please recalculate the cluster positions
1654 //Given the input event, loop over all the tracks, select the closest cluster as matched with fCutR
1655 //Store matched cluster indexes and residuals
1657 fMatchedTrackIndex ->Reset();
1658 fMatchedClusterIndex->Reset();
1659 fResidualPhi->Reset();
1660 fResidualEta->Reset();
1662 fMatchedTrackIndex ->Set(1000);
1663 fMatchedClusterIndex->Set(1000);
1664 fResidualPhi->Set(1000);
1665 fResidualEta->Set(1000);
1667 AliESDEvent* esdevent = dynamic_cast<AliESDEvent*> (event);
1668 AliAODEvent* aodevent = dynamic_cast<AliAODEvent*> (event);
1670 // init the magnetic field if not already on
1671 if(!TGeoGlobalMagField::Instance()->GetField())
1673 AliInfo("Init the magnetic field\n");
1676 esdevent->InitMagneticField();
1680 Double_t curSol = 30000*aodevent->GetMagneticField()/5.00668;
1681 Double_t curDip = 6000 *aodevent->GetMuonMagFieldScale();
1682 AliMagF *field = AliMagF::CreateFieldMap(curSol,curDip);
1683 TGeoGlobalMagField::Instance()->SetField(field);
1687 AliInfo("Mag Field not initialized, null esd/aod evetn pointers");
1692 TObjArray *clusterArray = 0x0;
1695 clusterArray = new TObjArray(event->GetNumberOfCaloClusters());
1696 for(Int_t icl=0; icl<event->GetNumberOfCaloClusters(); icl++)
1698 AliVCluster *cluster = (AliVCluster*) event->GetCaloCluster(icl);
1699 if(geom && !IsGoodCluster(cluster,geom,(AliVCaloCells*)event->GetEMCALCells())) continue;
1700 clusterArray->AddAt(cluster,icl);
1706 for (Int_t i=0; i<21;i++) cv[i]=0;
1707 for(Int_t itr=0; itr<event->GetNumberOfTracks(); itr++)
1709 AliExternalTrackParam *trackParam = 0;
1711 //If the input event is ESD, the starting point for extrapolation is TPCOut, if available, or TPCInner
1712 AliESDtrack *esdTrack = 0;
1713 AliAODTrack *aodTrack = 0;
1716 esdTrack = esdevent->GetTrack(itr);
1717 if(!esdTrack) continue;
1718 if(!IsAccepted(esdTrack)) continue;
1719 if(esdTrack->Pt()<fCutMinTrackPt) continue;
1720 Double_t phi = esdTrack->Phi()*TMath::RadToDeg();
1721 if(TMath::Abs(esdTrack->Eta())>0.8 || phi <= 20 || phi >= 240 ) continue;
1722 trackParam = const_cast<AliExternalTrackParam*>(esdTrack->GetInnerParam());
1725 //If the input event is AOD, the starting point for extrapolation is at vertex
1726 //AOD tracks are selected according to its filterbit.
1729 aodTrack = aodevent->GetTrack(itr);
1730 if(!aodTrack) continue;
1731 if(!aodTrack->TestFilterMask(fAODFilterMask)) continue; //Select AOD tracks that fulfill GetStandardITSTPCTrackCuts2010()
1732 if(aodTrack->Pt()<fCutMinTrackPt) continue;
1733 Double_t phi = aodTrack->Phi()*TMath::RadToDeg();
1734 if(TMath::Abs(aodTrack->Eta())>0.8 || phi <= 20 || phi >= 240 ) continue;
1735 Double_t pos[3],mom[3];
1736 aodTrack->GetXYZ(pos);
1737 aodTrack->GetPxPyPz(mom);
1738 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()));
1739 trackParam= new AliExternalTrackParam(pos,mom,cv,aodTrack->Charge());
1742 //Return if the input data is not "AOD" or "ESD"
1745 printf("Wrong input data type! Should be \"AOD\" or \"ESD\"\n");
1748 clusterArray->Clear();
1749 delete clusterArray;
1754 if(!trackParam) continue;
1756 //Extrapolate the track to EMCal surface
1757 AliExternalTrackParam emcalParam(*trackParam);
1759 if(!ExtrapolateTrackToEMCalSurface(&emcalParam, 430., fMass, fStepSurface, eta, phi))
1761 if(aodevent && trackParam) delete trackParam;
1767 // esdTrack->SetOuterParam(&emcalParam,AliExternalTrackParam::kMultSec);
1770 if(TMath::Abs(eta)>0.75 || (phi) < 70*TMath::DegToRad() || (phi) > 190*TMath::DegToRad())
1772 if(aodevent && trackParam) delete trackParam;
1777 //Find matched clusters
1779 Float_t dEta = -999, dPhi = -999;
1782 index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArray, dEta, dPhi);
1786 index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArr, dEta, dPhi);
1791 fMatchedTrackIndex ->AddAt(itr,matched);
1792 fMatchedClusterIndex ->AddAt(index,matched);
1793 fResidualEta ->AddAt(dEta,matched);
1794 fResidualPhi ->AddAt(dPhi,matched);
1797 if(aodevent && trackParam) delete trackParam;
1802 clusterArray->Clear();
1803 delete clusterArray;
1806 AliDebug(2,Form("Number of matched pairs = %d !\n",matched));
1808 fMatchedTrackIndex ->Set(matched);
1809 fMatchedClusterIndex ->Set(matched);
1810 fResidualPhi ->Set(matched);
1811 fResidualEta ->Set(matched);
1814 //________________________________________________________________________________
1815 Int_t AliEMCALRecoUtils::FindMatchedClusterInEvent(const AliESDtrack *track,
1816 const AliVEvent *event,
1817 const AliEMCALGeometry *geom,
1818 Float_t &dEta, Float_t &dPhi)
1821 // This function returns the index of matched cluster to input track
1822 // Returns -1 if no match is found
1824 Double_t phiV = track->Phi()*TMath::RadToDeg();
1825 if(TMath::Abs(track->Eta())>0.8 || phiV <= 20 || phiV >= 240 ) return index;
1826 AliExternalTrackParam *trackParam = const_cast<AliExternalTrackParam*>(track->GetInnerParam());
1827 if(!trackParam) return index;
1828 AliExternalTrackParam emcalParam(*trackParam);
1830 if(!ExtrapolateTrackToEMCalSurface(&emcalParam, 430., fMass, fStepSurface, eta, phi)) return index;
1831 if(TMath::Abs(eta)>0.75 || (phi) < 70*TMath::DegToRad() || (phi) > 190*TMath::DegToRad()) return index;
1833 TObjArray *clusterArr = new TObjArray(event->GetNumberOfCaloClusters());
1835 for(Int_t icl=0; icl<event->GetNumberOfCaloClusters(); icl++)
1837 AliVCluster *cluster = (AliVCluster*) event->GetCaloCluster(icl);
1838 if(geom && !IsGoodCluster(cluster,geom,(AliVCaloCells*)event->GetEMCALCells())) continue;
1839 clusterArr->AddAt(cluster,icl);
1842 index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArr, dEta, dPhi);
1843 clusterArr->Clear();
1849 //_______________________________________________________________________________________________
1850 Int_t AliEMCALRecoUtils::FindMatchedClusterInClusterArr(const AliExternalTrackParam *emcalParam,
1851 AliExternalTrackParam *trkParam,
1852 const TObjArray * clusterArr,
1853 Float_t &dEta, Float_t &dPhi)
1855 // Find matched cluster in array
1857 dEta=-999, dPhi=-999;
1858 Float_t dRMax = fCutR, dEtaMax=fCutEta, dPhiMax=fCutPhi;
1860 Float_t tmpEta=-999, tmpPhi=-999;
1862 Double_t exPos[3] = {0.,0.,0.};
1863 if(!emcalParam->GetXYZ(exPos)) return index;
1865 Float_t clsPos[3] = {0.,0.,0.};
1866 for(Int_t icl=0; icl<clusterArr->GetEntriesFast(); icl++)
1868 AliVCluster *cluster = dynamic_cast<AliVCluster*> (clusterArr->At(icl)) ;
1869 if(!cluster || !cluster->IsEMCAL()) continue;
1870 cluster->GetPosition(clsPos);
1871 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));
1872 if(dR > fClusterWindow) continue;
1874 AliExternalTrackParam trkPamTmp (*trkParam);//Retrieve the starting point every time before the extrapolation
1875 if(!ExtrapolateTrackToCluster(&trkPamTmp, cluster, fMass, fStepCluster, tmpEta, tmpPhi)) continue;
1878 Float_t tmpR=TMath::Sqrt(tmpEta*tmpEta + tmpPhi*tmpPhi);
1887 else if(fCutEtaPhiSeparate)
1889 if(TMath::Abs(tmpEta)<TMath::Abs(dEtaMax) && TMath::Abs(tmpPhi)<TMath::Abs(dPhiMax))
1898 printf("Error: please specify your cut criteria\n");
1899 printf("To cut on sqrt(dEta^2+dPhi^2), use: SwitchOnCutEtaPhiSum()\n");
1900 printf("To cut on dEta and dPhi separately, use: SwitchOnCutEtaPhiSeparate()\n");
1911 //------------------------------------------------------------------------------------
1912 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToEMCalSurface(AliExternalTrackParam *trkParam,
1913 const Double_t emcalR,
1914 const Double_t mass,
1915 const Double_t step,
1919 //Extrapolate track to EMCAL surface
1921 eta = -999, phi = -999;
1922 if(!trkParam) return kFALSE;
1923 if(!AliTrackerBase::PropagateTrackToBxByBz(trkParam, emcalR, mass, step, kTRUE, 0.8, -1)) return kFALSE;
1924 Double_t trkPos[3] = {0.,0.,0.};
1925 if(!trkParam->GetXYZ(trkPos)) return kFALSE;
1926 TVector3 trkPosVec(trkPos[0],trkPos[1],trkPos[2]);
1927 eta = trkPosVec.Eta();
1928 phi = trkPosVec.Phi();
1930 phi += 2*TMath::Pi();
1935 //-----------------------------------------------------------------------------------
1936 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToPosition(AliExternalTrackParam *trkParam,
1937 const Float_t *clsPos,
1944 //Return the residual by extrapolating a track param to a global position
1948 if(!trkParam) return kFALSE;
1949 Double_t trkPos[3] = {0.,0.,0.};
1950 TVector3 vec(clsPos[0],clsPos[1],clsPos[2]);
1951 Double_t alpha = ((int)(vec.Phi()*TMath::RadToDeg()/20)+0.5)*20*TMath::DegToRad();
1952 vec.RotateZ(-alpha); //Rotate the cluster to the local extrapolation coordinate system
1953 if(!AliTrackerBase::PropagateTrackToBxByBz(trkParam, vec.X(), mass, step,kTRUE, 0.8, -1)) return kFALSE;
1954 if(!trkParam->GetXYZ(trkPos)) return kFALSE; //Get the extrapolated global position
1956 TVector3 clsPosVec(clsPos[0],clsPos[1],clsPos[2]);
1957 TVector3 trkPosVec(trkPos[0],trkPos[1],trkPos[2]);
1959 // track cluster matching
1960 tmpPhi = clsPosVec.DeltaPhi(trkPosVec); // tmpPhi is between -pi and pi
1961 tmpEta = clsPosVec.Eta()-trkPosVec.Eta();
1966 //----------------------------------------------------------------------------------
1967 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToCluster(AliExternalTrackParam *trkParam,
1968 const AliVCluster *cluster,
1969 const Double_t mass,
1970 const Double_t step,
1975 //Return the residual by extrapolating a track param to a cluster
1979 if(!cluster || !trkParam) return kFALSE;
1981 Float_t clsPos[3] = {0.,0.,0.};
1982 cluster->GetPosition(clsPos);
1984 return ExtrapolateTrackToPosition(trkParam, clsPos, mass, step, tmpEta, tmpPhi);
1987 //---------------------------------------------------------------------------------
1988 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToCluster(AliExternalTrackParam *trkParam,
1989 const AliVCluster *cluster,
1994 //Return the residual by extrapolating a track param to a clusterfStepCluster
1997 return ExtrapolateTrackToCluster(trkParam, cluster, fMass, fStepCluster, tmpEta, tmpPhi);
2000 //_______________________________________________________________________
2001 void AliEMCALRecoUtils::GetMatchedResiduals(const Int_t clsIndex,
2002 Float_t &dEta, Float_t &dPhi)
2004 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2005 //Get the residuals dEta and dPhi for this cluster to the closest track
2006 //Works with ESDs and AODs
2008 if( FindMatchedPosForCluster(clsIndex) >= 999 )
2010 AliDebug(2,"No matched tracks found!\n");
2015 dEta = fResidualEta->At(FindMatchedPosForCluster(clsIndex));
2016 dPhi = fResidualPhi->At(FindMatchedPosForCluster(clsIndex));
2019 //______________________________________________________________________________________________
2020 void AliEMCALRecoUtils::GetMatchedClusterResiduals(Int_t trkIndex, Float_t &dEta, Float_t &dPhi)
2022 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2023 //Get the residuals dEta and dPhi for this track to the closest cluster
2024 //Works with ESDs and AODs
2026 if( FindMatchedPosForTrack(trkIndex) >= 999 )
2028 AliDebug(2,"No matched cluster found!\n");
2033 dEta = fResidualEta->At(FindMatchedPosForTrack(trkIndex));
2034 dPhi = fResidualPhi->At(FindMatchedPosForTrack(trkIndex));
2037 //__________________________________________________________
2038 Int_t AliEMCALRecoUtils::GetMatchedTrackIndex(Int_t clsIndex)
2040 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2041 //Get the index of matched track to this cluster
2042 //Works with ESDs and AODs
2044 if(IsClusterMatched(clsIndex))
2045 return fMatchedTrackIndex->At(FindMatchedPosForCluster(clsIndex));
2050 //__________________________________________________________
2051 Int_t AliEMCALRecoUtils::GetMatchedClusterIndex(Int_t trkIndex)
2053 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2054 //Get the index of matched cluster to this track
2055 //Works with ESDs and AODs
2057 if(IsTrackMatched(trkIndex))
2058 return fMatchedClusterIndex->At(FindMatchedPosForTrack(trkIndex));
2063 //______________________________________________________________
2064 Bool_t AliEMCALRecoUtils::IsClusterMatched(Int_t clsIndex) const
2066 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2067 //Returns if the cluster has a match
2068 if(FindMatchedPosForCluster(clsIndex) < 999)
2074 //____________________________________________________________
2075 Bool_t AliEMCALRecoUtils::IsTrackMatched(Int_t trkIndex) const
2077 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2078 //Returns if the track has a match
2079 if(FindMatchedPosForTrack(trkIndex) < 999)
2085 //______________________________________________________________________
2086 UInt_t AliEMCALRecoUtils::FindMatchedPosForCluster(Int_t clsIndex) const
2088 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2089 //Returns the position of the match in the fMatchedClusterIndex array
2090 Float_t tmpR = fCutR;
2093 for(Int_t i=0; i<fMatchedClusterIndex->GetSize(); i++)
2095 if(fMatchedClusterIndex->At(i)==clsIndex)
2097 Float_t r = TMath::Sqrt(fResidualEta->At(i)*fResidualEta->At(i) + fResidualPhi->At(i)*fResidualPhi->At(i));
2102 AliDebug(3,Form("Matched cluster index: index: %d, dEta: %2.4f, dPhi: %2.4f.\n",
2103 fMatchedClusterIndex->At(i),fResidualEta->At(i),fResidualPhi->At(i)));
2110 //____________________________________________________________________
2111 UInt_t AliEMCALRecoUtils::FindMatchedPosForTrack(Int_t trkIndex) const
2113 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2114 //Returns the position of the match in the fMatchedTrackIndex array
2115 Float_t tmpR = fCutR;
2118 for(Int_t i=0; i<fMatchedTrackIndex->GetSize(); i++)
2120 if(fMatchedTrackIndex->At(i)==trkIndex)
2122 Float_t r = TMath::Sqrt(fResidualEta->At(i)*fResidualEta->At(i) + fResidualPhi->At(i)*fResidualPhi->At(i));
2127 AliDebug(3,Form("Matched track index: index: %d, dEta: %2.4f, dPhi: %2.4f.\n",
2128 fMatchedTrackIndex->At(i),fResidualEta->At(i),fResidualPhi->At(i)));
2135 //__________________________________________________________________________
2136 Bool_t AliEMCALRecoUtils::IsGoodCluster(AliVCluster *cluster,
2137 const AliEMCALGeometry *geom,
2138 AliVCaloCells* cells,const Int_t bc)
2140 // check if the cluster survives some quality cut
2143 Bool_t isGood=kTRUE;
2145 if(!cluster || !cluster->IsEMCAL()) return kFALSE;
2147 if(ClusterContainsBadChannel(geom,cluster->GetCellsAbsId(),cluster->GetNCells())) return kFALSE;
2149 if(!CheckCellFiducialRegion(geom,cluster,cells)) return kFALSE;
2151 if(IsExoticCluster(cluster, cells,bc)) return kFALSE;
2156 //__________________________________________________________
2157 Bool_t AliEMCALRecoUtils::IsAccepted(AliESDtrack *esdTrack)
2159 // Given a esd track, return whether the track survive all the cuts
2161 // The different quality parameter are first
2162 // retrieved from the track. then it is found out what cuts the
2163 // track did not survive and finally the cuts are imposed.
2165 UInt_t status = esdTrack->GetStatus();
2167 Int_t nClustersITS = esdTrack->GetITSclusters(0);
2168 Int_t nClustersTPC = esdTrack->GetTPCclusters(0);
2170 Float_t chi2PerClusterITS = -1;
2171 Float_t chi2PerClusterTPC = -1;
2172 if (nClustersITS!=0)
2173 chi2PerClusterITS = esdTrack->GetITSchi2()/Float_t(nClustersITS);
2174 if (nClustersTPC!=0)
2175 chi2PerClusterTPC = esdTrack->GetTPCchi2()/Float_t(nClustersTPC);
2179 if(fTrackCutsType==kGlobalCut)
2181 Float_t maxDCAToVertexXYPtDep = 0.0182 + 0.0350/TMath::Power(esdTrack->Pt(),1.01); //This expression comes from AliESDtrackCuts::GetStandardITSTPCTrackCuts2010()
2182 //AliDebug(3,Form("Track pT = %f, DCAtoVertexXY = %f",esdTrack->Pt(),MaxDCAToVertexXYPtDep));
2183 SetMaxDCAToVertexXY(maxDCAToVertexXYPtDep); //Set pT dependent DCA cut to vertex in x-y plane
2189 esdTrack->GetImpactParameters(b,bCov);
2190 if (bCov[0]<=0 || bCov[2]<=0)
2192 AliDebug(1, "Estimated b resolution lower or equal zero!");
2193 bCov[0]=0; bCov[2]=0;
2196 Float_t dcaToVertexXY = b[0];
2197 Float_t dcaToVertexZ = b[1];
2198 Float_t dcaToVertex = -1;
2200 if (fCutDCAToVertex2D)
2201 dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY/fCutMaxDCAToVertexXY/fCutMaxDCAToVertexXY + dcaToVertexZ*dcaToVertexZ/fCutMaxDCAToVertexZ/fCutMaxDCAToVertexZ);
2203 dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY + dcaToVertexZ*dcaToVertexZ);
2207 Bool_t cuts[kNCuts];
2208 for (Int_t i=0; i<kNCuts; i++) cuts[i]=kFALSE;
2210 // track quality cuts
2211 if (fCutRequireTPCRefit && (status&AliESDtrack::kTPCrefit)==0)
2213 if (fCutRequireITSRefit && (status&AliESDtrack::kITSrefit)==0)
2215 if (nClustersTPC<fCutMinNClusterTPC)
2217 if (nClustersITS<fCutMinNClusterITS)
2219 if (chi2PerClusterTPC>fCutMaxChi2PerClusterTPC)
2221 if (chi2PerClusterITS>fCutMaxChi2PerClusterITS)
2223 if (!fCutAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0)
2225 if (fCutDCAToVertex2D && dcaToVertex > 1)
2227 if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexXY) > fCutMaxDCAToVertexXY)
2229 if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexZ) > fCutMaxDCAToVertexZ)
2232 if(fTrackCutsType==kGlobalCut)
2234 //Require at least one SPD point + anything else in ITS
2235 if( (esdTrack->HasPointOnITSLayer(0) || esdTrack->HasPointOnITSLayer(1)) == kFALSE)
2240 for (Int_t i=0; i<kNCuts; i++)
2241 if (cuts[i]) { cut = kTRUE ; }
2250 //_____________________________________
2251 void AliEMCALRecoUtils::InitTrackCuts()
2253 //Intilize the track cut criteria
2254 //By default these cuts are set according to AliESDtrackCuts::GetStandardTPCOnlyTrackCuts()
2255 //Also you can customize the cuts using the setters
2257 switch (fTrackCutsType)
2261 AliInfo(Form("Track cuts for matching: GetStandardTPCOnlyTrackCuts()"));
2263 SetMinNClustersTPC(70);
2264 SetMaxChi2PerClusterTPC(4);
2265 SetAcceptKinkDaughters(kFALSE);
2266 SetRequireTPCRefit(kFALSE);
2269 SetRequireITSRefit(kFALSE);
2270 SetMaxDCAToVertexZ(3.2);
2271 SetMaxDCAToVertexXY(2.4);
2272 SetDCAToVertex2D(kTRUE);
2279 AliInfo(Form("Track cuts for matching: GetStandardITSTPCTrackCuts2010(kTURE)"));
2281 SetMinNClustersTPC(70);
2282 SetMaxChi2PerClusterTPC(4);
2283 SetAcceptKinkDaughters(kFALSE);
2284 SetRequireTPCRefit(kTRUE);
2287 SetRequireITSRefit(kTRUE);
2288 SetMaxDCAToVertexZ(2);
2289 SetMaxDCAToVertexXY();
2290 SetDCAToVertex2D(kFALSE);
2297 AliInfo(Form("Track cuts for matching: Loose cut w/o DCA cut"));
2298 SetMinNClustersTPC(50);
2299 SetAcceptKinkDaughters(kTRUE);
2307 //________________________________________________________________________
2308 void AliEMCALRecoUtils::SetClusterMatchedToTrack(const AliESDEvent *event)
2310 // Checks if tracks are matched to EMC clusters and set the matched EMCAL cluster index to ESD track.
2312 Int_t nTracks = event->GetNumberOfTracks();
2313 for (Int_t iTrack = 0; iTrack < nTracks; ++iTrack)
2315 AliESDtrack* track = event->GetTrack(iTrack);
2318 AliWarning(Form("Could not receive track %d", iTrack));
2322 Int_t matchClusIndex = GetMatchedClusterIndex(iTrack);
2323 track->SetEMCALcluster(matchClusIndex); //sets -1 if track not matched within residual
2324 if(matchClusIndex != -1)
2325 track->SetStatus(AliESDtrack::kEMCALmatch);
2327 track->ResetStatus(AliESDtrack::kEMCALmatch);
2329 AliDebug(2,"Track matched to closest cluster");
2332 //_________________________________________________________________________
2333 void AliEMCALRecoUtils::SetTracksMatchedToCluster(const AliESDEvent *event)
2335 // Checks if EMC clusters are matched to ESD track.
2336 // Adds track indexes of all the tracks matched to a cluster withing residuals in ESDCalocluster.
2338 for (Int_t iClus=0; iClus < event->GetNumberOfCaloClusters(); ++iClus)
2340 AliESDCaloCluster *cluster = event->GetCaloCluster(iClus);
2341 if (!cluster->IsEMCAL())
2344 Int_t nTracks = event->GetNumberOfTracks();
2345 TArrayI arrayTrackMatched(nTracks);
2347 // Get the closest track matched to the cluster
2349 Int_t matchTrackIndex = GetMatchedTrackIndex(iClus);
2350 if (matchTrackIndex != -1)
2352 arrayTrackMatched[nMatched] = matchTrackIndex;
2356 // Get all other tracks matched to the cluster
2357 for(Int_t iTrk=0; iTrk<nTracks; ++iTrk)
2359 AliESDtrack* track = event->GetTrack(iTrk);
2360 if(iTrk == matchTrackIndex) continue;
2361 if(track->GetEMCALcluster() == iClus)
2363 arrayTrackMatched[nMatched] = iTrk;
2368 //printf("Tender::SetTracksMatchedToCluster - cluster E %f, N matches %d, first match %d\n",cluster->E(),nMatched,arrayTrackMatched[0]);
2370 arrayTrackMatched.Set(nMatched);
2371 cluster->AddTracksMatched(arrayTrackMatched);
2373 Float_t eta= -999, phi = -999;
2374 if (matchTrackIndex != -1)
2375 GetMatchedResiduals(iClus, eta, phi);
2376 cluster->SetTrackDistance(phi, eta);
2379 AliDebug(2,"Cluster matched to tracks");
2383 //___________________________________________________
2384 void AliEMCALRecoUtils::Print(const Option_t *) const
2388 printf("AliEMCALRecoUtils Settings: \n");
2389 printf("Misalignment shifts\n");
2390 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,
2391 fMisalTransShift[i*3],fMisalTransShift[i*3+1],fMisalTransShift[i*3+2],
2392 fMisalRotShift[i*3], fMisalRotShift[i*3+1], fMisalRotShift[i*3+2] );
2393 printf("Non linearity function %d, parameters:\n", fNonLinearityFunction);
2394 for(Int_t i=0; i<6; i++) printf("param[%d]=%f\n",i, fNonLinearityParams[i]);
2396 printf("Position Recalculation option %d, Particle Type %d, fW0 %2.2f, Recalibrate Data %d \n",fPosAlgo,fParticleType,fW0, fRecalibration);
2398 printf("Matching criteria: ");
2401 printf("sqrt(dEta^2+dPhi^2)<%4.3f\n",fCutR);
2403 else if(fCutEtaPhiSeparate)
2405 printf("dEta<%4.3f, dPhi<%4.3f\n",fCutEta,fCutPhi);
2410 printf("please specify your cut criteria\n");
2411 printf("To cut on sqrt(dEta^2+dPhi^2), use: SwitchOnCutEtaPhiSum()\n");
2412 printf("To cut on dEta and dPhi separately, use: SwitchOnCutEtaPhiSeparate()\n");
2415 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);
2416 printf("Cluster selection window: dR < %2.0f\n",fClusterWindow);
2418 printf("Track cuts: \n");
2419 printf("Minimum track pT: %1.2f\n",fCutMinTrackPt);
2420 printf("AOD track selection mask: %d\n",fAODFilterMask);
2421 printf("TPCRefit = %d, ITSRefit = %d\n",fCutRequireTPCRefit,fCutRequireITSRefit);
2422 printf("AcceptKinks = %d\n",fCutAcceptKinkDaughters);
2423 printf("MinNCulsterTPC = %d, MinNClusterITS = %d\n",fCutMinNClusterTPC,fCutMinNClusterITS);
2424 printf("MaxChi2TPC = %2.2f, MaxChi2ITS = %2.2f\n",fCutMaxChi2PerClusterTPC,fCutMaxChi2PerClusterITS);
2425 printf("DCSToVertex2D = %d, MaxDCAToVertexXY = %2.2f, MaxDCAToVertexZ = %2.2f\n",fCutDCAToVertex2D,fCutMaxDCAToVertexXY,fCutMaxDCAToVertexZ);
2428 //_____________________________________________________________________
2429 void AliEMCALRecoUtils::SetRunDependentCorrections(Int_t runnumber)
2431 //Get EMCAL time dependent corrections from file and put them in the recalibration histograms
2432 //Do it only once and only if it is requested
2434 if(!fUseRunCorrectionFactors) return;
2435 if(fRunCorrectionFactorsSet) return;
2437 AliInfo(Form("AliEMCALRecoUtils::GetRunDependentCorrections() - Get Correction Factors for Run number %d\n",runnumber));
2439 AliEMCALCalibTimeDepCorrection *corr = new AliEMCALCalibTimeDepCorrection();
2440 corr->ReadRootInfo(Form("CorrectionFiles/Run%d_Correction.root",runnumber));
2442 SwitchOnRecalibration();
2443 for(Int_t ism = 0; ism < 4; ism++)
2445 for(Int_t icol = 0; icol < 48; icol++)
2447 for(Int_t irow = 0; irow < 24; irow++)
2449 Float_t orgRecalFactor = GetEMCALChannelRecalibrationFactors(ism)->GetBinContent(icol,irow);
2450 Float_t newRecalFactor = orgRecalFactor*corr->GetCorrection(ism, icol,irow,0);
2451 GetEMCALChannelRecalibrationFactors(ism)->SetBinContent(icol,irow,newRecalFactor);
2452 //printf("ism %d, icol %d, irow %d, corrections : org %f, time dep %f, final %f (org*time %f)\n",ism, icol, irow,
2453 // orgRecalFactor, corr->GetCorrection(ism, icol,irow,0),
2454 // (GetEMCALChannelRecalibrationFactors(ism))->GetBinContent(icol,irow),newRecalFactor);
2458 fRunCorrectionFactorsSet = kTRUE;