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){
221 // assign or copy construct
223 *fResidualEta = *reco.fResidualEta;
225 else fResidualEta = new TArrayF(*reco.fResidualEta);
228 if(fResidualEta)delete fResidualEta;
232 if(reco.fResidualPhi){
233 // assign or copy construct
235 *fResidualPhi = *reco.fResidualPhi;
237 else fResidualPhi = new TArrayF(*reco.fResidualPhi);
240 if(fResidualPhi)delete fResidualPhi;
244 if(reco.fMatchedTrackIndex){
245 // assign or copy construct
246 if(fMatchedTrackIndex){
247 *fMatchedTrackIndex = *reco.fMatchedTrackIndex;
249 else fMatchedTrackIndex = new TArrayI(*reco.fMatchedTrackIndex);
252 if(fMatchedTrackIndex)delete fMatchedTrackIndex;
253 fMatchedTrackIndex = 0;
256 if(reco.fMatchedClusterIndex){
257 // assign or copy construct
258 if(fMatchedClusterIndex){
259 *fMatchedClusterIndex = *reco.fMatchedClusterIndex;
261 else fMatchedClusterIndex = new TArrayI(*reco.fMatchedClusterIndex);
264 if(fMatchedClusterIndex)delete fMatchedClusterIndex;
265 fMatchedClusterIndex = 0;
272 //_____________________________________
273 AliEMCALRecoUtils::~AliEMCALRecoUtils()
277 if(fEMCALRecalibrationFactors) {
278 fEMCALRecalibrationFactors->Clear();
279 delete fEMCALRecalibrationFactors;
282 if(fEMCALTimeRecalibrationFactors) {
283 fEMCALTimeRecalibrationFactors->Clear();
284 delete fEMCALTimeRecalibrationFactors;
287 if(fEMCALBadChannelMap) {
288 fEMCALBadChannelMap->Clear();
289 delete fEMCALBadChannelMap;
292 delete fMatchedTrackIndex ;
293 delete fMatchedClusterIndex ;
294 delete fResidualEta ;
295 delete fResidualPhi ;
301 //_______________________________________________________________________________
302 Bool_t AliEMCALRecoUtils::AcceptCalibrateCell(const Int_t absID, const Int_t bc,
303 Float_t & amp, Double_t & time,
304 AliVCaloCells* cells)
306 // Reject cell if criteria not passed and calibrate it
308 AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
310 if(absID < 0 || absID >= 24*48*geom->GetNumberOfSuperModules()) return kFALSE;
312 Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1;
314 if(!geom->GetCellIndex(absID,imod,iTower,iIphi,iIeta))
316 // cell absID does not exist
321 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);
323 // Do not include bad channels found in analysis,
324 if( IsBadChannelsRemovalSwitchedOn() && GetEMCALChannelStatus(imod, ieta, iphi))
330 amp = cells->GetCellAmplitude(absID);
331 if(!fCellsRecalibrated && IsRecalibrationOn())
332 amp *= GetEMCALChannelRecalibrationFactor(imod,ieta,iphi);
336 time = cells->GetCellTime(absID);
338 RecalibrateCellTime(absID,bc,time);
343 //_____________________________________________________________________________
344 Bool_t AliEMCALRecoUtils::CheckCellFiducialRegion(const AliEMCALGeometry* geom,
345 const AliVCluster* cluster,
346 AliVCaloCells* cells)
348 // Given the list of AbsId of the cluster, get the maximum cell and
349 // check if there are fNCellsFromBorder from the calorimeter border
352 AliInfo("Cluster pointer null!");
356 //If the distance to the border is 0 or negative just exit accept all clusters
357 if(cells->GetType()==AliVCaloCells::kEMCALCell && fNCellsFromEMCALBorder <= 0 ) return kTRUE;
359 Int_t absIdMax = -1, iSM =-1, ieta = -1, iphi = -1;
360 Bool_t shared = kFALSE;
361 GetMaxEnergyCell(geom, cells, cluster, absIdMax, iSM, ieta, iphi, shared);
363 AliDebug(2,Form("Cluster Max AbsId %d, Cell Energy %2.2f, Cluster Energy %2.2f, Ncells from border %d, EMCAL eta=0 %d\n",
364 absIdMax, cells->GetCellAmplitude(absIdMax), cluster->E(), fNCellsFromEMCALBorder, fNoEMCALBorderAtEta0));
366 if(absIdMax==-1) return kFALSE;
368 //Check if the cell is close to the borders:
369 Bool_t okrow = kFALSE;
370 Bool_t okcol = kFALSE;
372 if(iSM < 0 || iphi < 0 || ieta < 0 ) {
373 AliFatal(Form("Negative value for super module: %d, or cell ieta: %d, or cell iphi: %d, check EMCAL geometry name\n",
379 if(iphi >= fNCellsFromEMCALBorder && iphi < 24-fNCellsFromEMCALBorder) okrow =kTRUE;
381 else if (iSM >=10 && ( ( geom->GetEMCGeometry()->GetGeoName()).Contains("12SMV1"))) {
382 if(iphi >= fNCellsFromEMCALBorder && iphi < 8-fNCellsFromEMCALBorder) okrow =kTRUE; //1/3 sm case
385 if(iphi >= fNCellsFromEMCALBorder && iphi < 12-fNCellsFromEMCALBorder) okrow =kTRUE; // half SM case
389 if(!fNoEMCALBorderAtEta0){
390 if(ieta > fNCellsFromEMCALBorder && ieta < 48-fNCellsFromEMCALBorder) okcol =kTRUE;
394 if(ieta >= fNCellsFromEMCALBorder) okcol = kTRUE;
397 if(ieta < 48-fNCellsFromEMCALBorder) okcol = kTRUE;
401 AliDebug(2,Form("EMCAL Cluster in %d cells fiducial volume: ieta %d, iphi %d, SM %d: column? %d, row? %d\nq",
402 fNCellsFromEMCALBorder, ieta, iphi, iSM, okcol, okrow));
404 if (okcol && okrow) {
405 //printf("Accept\n");
409 //printf("Reject\n");
410 AliDebug(2,Form("Reject cluster in border, max cell : ieta %d, iphi %d, SM %d\n",ieta, iphi, iSM));
417 //_______________________________________________________________________________
418 Bool_t AliEMCALRecoUtils::ClusterContainsBadChannel(const AliEMCALGeometry* geom,
419 const UShort_t* cellList,
422 // Check that in the cluster cells, there is no bad channel of those stored
423 // in fEMCALBadChannelMap or fPHOSBadChannelMap
425 if(!fRemoveBadChannels) return kFALSE;
426 if(!fEMCALBadChannelMap) return kFALSE;
431 for(Int_t iCell = 0; iCell<nCells; iCell++){
433 //Get the column and row
434 Int_t iTower = -1, iIphi = -1, iIeta = -1;
435 geom->GetCellIndex(cellList[iCell],imod,iTower,iIphi,iIeta);
436 if(fEMCALBadChannelMap->GetEntries() <= imod) continue;
437 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,irow,icol);
438 if(GetEMCALChannelStatus(imod, icol, irow)){
439 AliDebug(2,Form("Cluster with bad channel: SM %d, col %d, row %d\n",imod, icol, irow));
443 }// cell cluster loop
448 //_____________________________________________________________________________________________
449 Bool_t AliEMCALRecoUtils::IsExoticCell(const Int_t absID, AliVCaloCells* cells, const Int_t bc)
451 // Look to cell neighbourhood and reject if it seems exotic
452 // Do before recalibrating the cells
454 if(!fRejectExoticCells) return kFALSE;
456 AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
458 Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1;
459 geom->GetCellIndex(absID,imod,iTower,iIphi,iIeta);
460 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);
462 //Get close cells index, energy and time, not in corners
463 Int_t absID1 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi+1, ieta);
464 Int_t absID2 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi-1, ieta);
465 Int_t absID3 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta+1);
466 Int_t absID4 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta-1);
468 Float_t ecell = 0, ecell1 = 0, ecell2 = 0, ecell3 = 0, ecell4 = 0;
469 Double_t tcell = 0, tcell1 = 0, tcell2 = 0, tcell3 = 0, tcell4 = 0;
470 Bool_t accept = 0, accept1 = 0, accept2 = 0, accept3 = 0, accept4 = 0;
472 accept = AcceptCalibrateCell(absID, bc, ecell ,tcell ,cells);
474 if(!accept) return kTRUE; // reject this cell
476 if(ecell < fExoticCellMinAmplitude) return kFALSE; // do not reject low energy cells
478 accept1 = AcceptCalibrateCell(absID1,bc, ecell1,tcell1,cells);
479 accept2 = AcceptCalibrateCell(absID2,bc, ecell2,tcell2,cells);
480 accept3 = AcceptCalibrateCell(absID3,bc, ecell3,tcell3,cells);
481 accept4 = AcceptCalibrateCell(absID4,bc, ecell4,tcell4,cells);
484 printf("Cell absID %d \n",absID);
485 printf("\t accept1 %d, accept2 %d, accept3 %d, accept4 %d\n",
486 accept1,accept2,accept3,accept4);
487 printf("\t id %d: id1 %d, id2 %d, id3 %d, id4 %d\n",
488 absID,absID1,absID2,absID3,absID4);
489 printf("\t e %f: e1 %f, e2 %f, e3 %f, e4 %f\n",
490 ecell,ecell1,ecell2,ecell3,ecell4);
491 printf("\t t %f: t1 %f, t2 %f, t3 %f, t4 %f;\n dt1 %f, dt2 %f, dt3 %f, dt4 %f\n",
492 tcell*1.e9,tcell1*1.e9,tcell2*1.e9,tcell3*1.e9,tcell4*1.e9,
493 TMath::Abs(tcell-tcell1)*1.e9, TMath::Abs(tcell-tcell2)*1.e9, TMath::Abs(tcell-tcell3)*1.e9, TMath::Abs(tcell-tcell4)*1.e9);
496 if(TMath::Abs(tcell-tcell1)*1.e9 > fExoticCellDiffTime) ecell1 = 0 ;
497 if(TMath::Abs(tcell-tcell2)*1.e9 > fExoticCellDiffTime) ecell2 = 0 ;
498 if(TMath::Abs(tcell-tcell3)*1.e9 > fExoticCellDiffTime) ecell3 = 0 ;
499 if(TMath::Abs(tcell-tcell4)*1.e9 > fExoticCellDiffTime) ecell4 = 0 ;
501 Float_t eCross = ecell1+ecell2+ecell3+ecell4;
503 //printf("\t eCell %f, eCross %f, 1-eCross/eCell %f\n",ecell,eCross,1-eCross/ecell);
505 if(1-eCross/ecell > fExoticCellFraction) {
506 AliDebug(2,Form("AliEMCALRecoUtils::IsExoticCell() - EXOTIC CELL id %d, eCell %f, eCross %f, 1-eCross/eCell %f\n",
507 absID,ecell,eCross,1-eCross/ecell));
514 //___________________________________________________________________
515 Bool_t AliEMCALRecoUtils::IsExoticCluster(const AliVCluster *cluster,
516 AliVCaloCells *cells,
519 // Check if the cluster highest energy tower is exotic
522 AliInfo("Cluster pointer null!");
526 if(!fRejectExoticCluster) return kFALSE;
528 // Get highest energy tower
529 AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
530 Int_t iSupMod = -1, absId = -1, ieta = -1, iphi = -1;
531 Bool_t shared = kFALSE;
532 GetMaxEnergyCell(geom, cells, cluster, absId, iSupMod, ieta, iphi, shared);
534 return IsExoticCell(absId,cells,bc);
537 //_______________________________________________________________________
538 Float_t AliEMCALRecoUtils::SmearClusterEnergy(const AliVCluster* cluster)
540 //In case of MC analysis, smear energy to match resolution/calibration in real data
543 AliInfo("Cluster pointer null!");
547 Float_t energy = cluster->E() ;
548 Float_t rdmEnergy = energy ;
549 if(fSmearClusterEnergy){
550 rdmEnergy = fRandom.Gaus(energy,fSmearClusterParam[0] * TMath::Sqrt(energy) +
551 fSmearClusterParam[1] * energy +
552 fSmearClusterParam[2] );
553 AliDebug(2, Form("Energy: original %f, smeared %f\n", energy, rdmEnergy));
559 //____________________________________________________________________________
560 Float_t AliEMCALRecoUtils::CorrectClusterEnergyLinearity(AliVCluster* cluster)
562 // Correct cluster energy from non linearity functions
565 AliInfo("Cluster pointer null!");
569 Float_t energy = cluster->E();
571 switch (fNonLinearityFunction) {
575 //Non-Linearity correction (from MC with function ([0]*exp(-[1]/E))+(([2]/([3]*2.*TMath::Pi())*exp(-(E-[4])^2/(2.*[3]^2)))))
576 //Double_t fNonLinearityParams[0] = 1.014;
577 //Double_t fNonLinearityParams[1] = -0.03329;
578 //Double_t fNonLinearityParams[2] = -0.3853;
579 //Double_t fNonLinearityParams[3] = 0.5423;
580 //Double_t fNonLinearityParams[4] = -0.4335;
581 energy *= (fNonLinearityParams[0]*exp(-fNonLinearityParams[1]/energy))+
582 ((fNonLinearityParams[2]/(fNonLinearityParams[3]*2.*TMath::Pi())*
583 exp(-(energy-fNonLinearityParams[4])*(energy-fNonLinearityParams[4])/(2.*fNonLinearityParams[3]*fNonLinearityParams[3]))));
589 //Non-Linearity correction (from Olga Data with function p0+p1*exp(-p2*E))
590 //Double_t fNonLinearityParams[0] = 1.04;
591 //Double_t fNonLinearityParams[1] = -0.1445;
592 //Double_t fNonLinearityParams[2] = 1.046;
593 energy /= (fNonLinearityParams[0]+fNonLinearityParams[1]*exp(-fNonLinearityParams[2]*energy)); //Olga function
597 case kPi0GammaConversion:
599 //Non-Linearity correction (Nicolas from Dimitri Data with function C*[1-a*exp(-b*E)])
600 //fNonLinearityParams[0] = 0.139393/0.1349766;
601 //fNonLinearityParams[1] = 0.0566186;
602 //fNonLinearityParams[2] = 0.982133;
603 energy /= fNonLinearityParams[0]*(1-fNonLinearityParams[1]*exp(-fNonLinearityParams[2]*energy));
610 //From beam test, Alexei's results, for different ZS thresholds
611 // th=30 MeV; th = 45 MeV; th = 75 MeV
612 //fNonLinearityParams[0] = 1.007; 1.003; 1.002
613 //fNonLinearityParams[1] = 0.894; 0.719; 0.797
614 //fNonLinearityParams[2] = 0.246; 0.334; 0.358
615 //Rescale the param[0] with 1.03
616 energy /= fNonLinearityParams[0]/(1+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]));
621 case kBeamTestCorrected:
623 //From beam test, corrected for material between beam and EMCAL
624 //fNonLinearityParams[0] = 0.99078
625 //fNonLinearityParams[1] = 0.161499;
626 //fNonLinearityParams[2] = 0.655166;
627 //fNonLinearityParams[3] = 0.134101;
628 //fNonLinearityParams[4] = 163.282;
629 //fNonLinearityParams[5] = 23.6904;
630 //fNonLinearityParams[6] = 0.978;
631 energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
637 AliDebug(2,"No correction on the energy\n");
645 //__________________________________________________
646 void AliEMCALRecoUtils::InitNonLinearityParam()
648 //Initialising Non Linearity Parameters
650 if(fNonLinearityFunction == kPi0MC) {
651 fNonLinearityParams[0] = 1.014;
652 fNonLinearityParams[1] = -0.03329;
653 fNonLinearityParams[2] = -0.3853;
654 fNonLinearityParams[3] = 0.5423;
655 fNonLinearityParams[4] = -0.4335;
658 if(fNonLinearityFunction == kPi0GammaGamma) {
659 fNonLinearityParams[0] = 1.04;
660 fNonLinearityParams[1] = -0.1445;
661 fNonLinearityParams[2] = 1.046;
664 if(fNonLinearityFunction == kPi0GammaConversion) {
665 fNonLinearityParams[0] = 0.139393;
666 fNonLinearityParams[1] = 0.0566186;
667 fNonLinearityParams[2] = 0.982133;
670 if(fNonLinearityFunction == kBeamTest) {
671 if(fNonLinearThreshold == 30) {
672 fNonLinearityParams[0] = 1.007;
673 fNonLinearityParams[1] = 0.894;
674 fNonLinearityParams[2] = 0.246;
676 if(fNonLinearThreshold == 45) {
677 fNonLinearityParams[0] = 1.003;
678 fNonLinearityParams[1] = 0.719;
679 fNonLinearityParams[2] = 0.334;
681 if(fNonLinearThreshold == 75) {
682 fNonLinearityParams[0] = 1.002;
683 fNonLinearityParams[1] = 0.797;
684 fNonLinearityParams[2] = 0.358;
688 if(fNonLinearityFunction == kBeamTestCorrected) {
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;
699 //_________________________________________________________
700 Float_t AliEMCALRecoUtils::GetDepth(const Float_t energy,
701 const Int_t iParticle,
702 const Int_t iSM) const
704 //Calculate shower depth for a given cluster energy and particle type
714 depth = x0 * (TMath::Log(energy*1000/ ecr) + 0.5); //Multiply energy by 1000 to transform to MeV
718 depth = x0 * (TMath::Log(energy*1000/ ecr) - 0.5); //Multiply energy by 1000 to transform to MeV
725 gGeoManager->cd("ALIC_1/XEN1_1");
726 TGeoNode *geoXEn1 = gGeoManager->GetCurrentNode();
727 TGeoNodeMatrix *geoSM = dynamic_cast<TGeoNodeMatrix *>(geoXEn1->GetDaughter(iSM));
729 TGeoVolume *geoSMVol = geoSM->GetVolume();
730 TGeoShape *geoSMShape = geoSMVol->GetShape();
731 TGeoBBox *geoBox = dynamic_cast<TGeoBBox *>(geoSMShape);
732 if(geoBox) depth = 0.5 * geoBox->GetDX()*2 ;
733 else AliFatal("Null GEANT box");
734 }else AliFatal("NULL GEANT node matrix");
737 depth = x0 * (TMath::Log(energy*1000 / ecr) - 0.5); //Multiply energy by 1000 to transform to MeV
743 depth = x0 * (TMath::Log(energy*1000 / ecr) + 0.5); //Multiply energy by 1000 to transform to MeV
749 //____________________________________________________________________
750 void AliEMCALRecoUtils::GetMaxEnergyCell(const AliEMCALGeometry *geom,
751 AliVCaloCells* cells,
752 const AliVCluster* clu,
759 //For a given CaloCluster gets the absId of the cell
760 //with maximum energy deposit.
763 Double_t eCell = -1.;
764 Float_t fraction = 1.;
765 Float_t recalFactor = 1.;
766 Int_t cellAbsId = -1 ;
774 AliInfo("Cluster pointer null!");
775 absId=-1; iSupMod0=-1, ieta = -1; iphi = -1; shared = -1;
779 for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++) {
780 cellAbsId = clu->GetCellAbsId(iDig);
781 fraction = clu->GetCellAmplitudeFraction(iDig);
782 //printf("a Cell %d, id, %d, amp %f, fraction %f\n",iDig,cellAbsId,cells->GetCellAmplitude(cellAbsId),fraction);
783 if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
784 geom->GetCellIndex(cellAbsId,iSupMod,iTower,iIphi,iIeta);
785 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);
786 if(iDig==0) iSupMod0=iSupMod;
787 else if(iSupMod0!=iSupMod) {
789 //printf("AliEMCALRecoUtils::GetMaxEnergyCell() - SHARED CLUSTER\n");
791 if(!fCellsRecalibrated && IsRecalibrationOn()) {
792 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
794 eCell = cells->GetCellAmplitude(cellAbsId)*fraction*recalFactor;
795 //printf("b Cell %d, id, %d, amp %f, fraction %f\n",iDig,cellAbsId,eCell,fraction);
799 //printf("\t new max: cell %d, e %f, ecell %f\n",maxId, eMax,eCell);
803 //Get from the absid the supermodule, tower and eta/phi numbers
804 geom->GetCellIndex(absId,iSupMod,iTower,iIphi,iIeta);
805 //Gives SuperModule and Tower numbers
806 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,
807 iIphi, iIeta,iphi,ieta);
808 //printf("Max id %d, iSM %d, col %d, row %d\n",absId,iSupMod,ieta,iphi);
809 //printf("Max end---\n");
812 //______________________________________
813 void AliEMCALRecoUtils::InitParameters()
815 // Initialize data members with default values
817 fParticleType = kPhoton;
818 fPosAlgo = kUnchanged;
821 fNonLinearityFunction = kNoCorrection;
822 fNonLinearThreshold = 30;
824 fExoticCellFraction = 0.97;
825 fExoticCellDiffTime = 1e6;
826 fExoticCellMinAmplitude = 0.5;
830 fCutEtaPhiSum = kTRUE;
831 fCutEtaPhiSeparate = kFALSE;
837 fClusterWindow = 100;
842 fTrackCutsType = kLooseCut;
845 fCutMinNClusterTPC = -1;
846 fCutMinNClusterITS = -1;
848 fCutMaxChi2PerClusterTPC = 1e10;
849 fCutMaxChi2PerClusterITS = 1e10;
851 fCutRequireTPCRefit = kFALSE;
852 fCutRequireITSRefit = kFALSE;
853 fCutAcceptKinkDaughters = kFALSE;
855 fCutMaxDCAToVertexXY = 1e10;
856 fCutMaxDCAToVertexZ = 1e10;
857 fCutDCAToVertex2D = kFALSE;
860 //Misalignment matrices
861 for(Int_t i = 0; i < 15 ; i++) {
862 fMisalTransShift[i] = 0.;
863 fMisalRotShift[i] = 0.;
867 for(Int_t i = 0; i < 7 ; i++) fNonLinearityParams[i] = 0.;
869 //For kBeamTestCorrected case, but default is no correction
870 fNonLinearityParams[0] = 0.99078;
871 fNonLinearityParams[1] = 0.161499;
872 fNonLinearityParams[2] = 0.655166;
873 fNonLinearityParams[3] = 0.134101;
874 fNonLinearityParams[4] = 163.282;
875 fNonLinearityParams[5] = 23.6904;
876 fNonLinearityParams[6] = 0.978;
878 //For kPi0GammaGamma case
879 //fNonLinearityParams[0] = 0.1457/0.1349766/1.038;
880 //fNonLinearityParams[1] = -0.02024/0.1349766/1.038;
881 //fNonLinearityParams[2] = 1.046;
883 //Cluster energy smearing
884 fSmearClusterEnergy = kFALSE;
885 fSmearClusterParam[0] = 0.07; // * sqrt E term
886 fSmearClusterParam[1] = 0.00; // * E term
887 fSmearClusterParam[2] = 0.00; // constant
890 //_____________________________________________________
891 void AliEMCALRecoUtils::InitEMCALRecalibrationFactors()
893 //Init EMCAL recalibration factors
894 AliDebug(2,"AliCalorimeterUtils::InitEMCALRecalibrationFactors()");
895 //In order to avoid rewriting the same histograms
896 Bool_t oldStatus = TH1::AddDirectoryStatus();
897 TH1::AddDirectory(kFALSE);
899 fEMCALRecalibrationFactors = new TObjArray(12);
900 for (int i = 0; i < 12; i++)
901 fEMCALRecalibrationFactors->Add(new TH2F(Form("EMCALRecalFactors_SM%d",i),
902 Form("EMCALRecalFactors_SM%d",i), 48, 0, 48, 24, 0, 24));
903 //Init the histograms with 1
904 for (Int_t sm = 0; sm < 12; sm++) {
905 for (Int_t i = 0; i < 48; i++) {
906 for (Int_t j = 0; j < 24; j++) {
907 SetEMCALChannelRecalibrationFactor(sm,i,j,1.);
911 fEMCALRecalibrationFactors->SetOwner(kTRUE);
912 fEMCALRecalibrationFactors->Compress();
914 //In order to avoid rewriting the same histograms
915 TH1::AddDirectory(oldStatus);
918 //_________________________________________________________
919 void AliEMCALRecoUtils::InitEMCALTimeRecalibrationFactors()
921 //Init EMCAL recalibration factors
922 AliDebug(2,"AliCalorimeterUtils::InitEMCALRecalibrationFactors()");
923 //In order to avoid rewriting the same histograms
924 Bool_t oldStatus = TH1::AddDirectoryStatus();
925 TH1::AddDirectory(kFALSE);
927 fEMCALTimeRecalibrationFactors = new TObjArray(4);
928 for (int i = 0; i < 4; i++)
929 fEMCALTimeRecalibrationFactors->Add(new TH1F(Form("hAllTimeAvBC%d",i),
930 Form("hAllTimeAvBC%d",i),
931 48*24*12,0.,48*24*12) );
932 //Init the histograms with 1
933 for (Int_t bc = 0; bc < 4; bc++) {
934 for (Int_t i = 0; i < 48*24*12; i++)
935 SetEMCALChannelTimeRecalibrationFactor(bc,i,0.);
938 fEMCALTimeRecalibrationFactors->SetOwner(kTRUE);
939 fEMCALTimeRecalibrationFactors->Compress();
941 //In order to avoid rewriting the same histograms
942 TH1::AddDirectory(oldStatus);
945 //____________________________________________________
946 void AliEMCALRecoUtils::InitEMCALBadChannelStatusMap()
948 //Init EMCAL bad channels map
949 AliDebug(2,"AliEMCALRecoUtils::InitEMCALBadChannelStatusMap()");
950 //In order to avoid rewriting the same histograms
951 Bool_t oldStatus = TH1::AddDirectoryStatus();
952 TH1::AddDirectory(kFALSE);
954 fEMCALBadChannelMap = new TObjArray(12);
955 //TH2F * hTemp = new TH2I("EMCALBadChannelMap","EMCAL SuperModule bad channel map", 48, 0, 48, 24, 0, 24);
956 for (int i = 0; i < 12; i++) {
957 fEMCALBadChannelMap->Add(new TH2I(Form("EMCALBadChannelMap_Mod%d",i),Form("EMCALBadChannelMap_Mod%d",i), 48, 0, 48, 24, 0, 24));
960 fEMCALBadChannelMap->SetOwner(kTRUE);
961 fEMCALBadChannelMap->Compress();
963 //In order to avoid rewriting the same histograms
964 TH1::AddDirectory(oldStatus);
967 //____________________________________________________________________________
968 void AliEMCALRecoUtils::RecalibrateClusterEnergy(const AliEMCALGeometry* geom,
969 AliVCluster * cluster,
970 AliVCaloCells * cells,
973 // Recalibrate the cluster energy and Time, considering the recalibration map
974 // and the energy of the cells and time that compose the cluster.
975 // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
978 AliInfo("Cluster pointer null!");
982 //Get the cluster number of cells and list of absId, check what kind of cluster do we have.
983 UShort_t * index = cluster->GetCellsAbsId() ;
984 Double_t * fraction = cluster->GetCellsAmplitudeFraction() ;
985 Int_t ncells = cluster->GetNCells();
987 //Initialize some used variables
990 Int_t icol =-1, irow =-1, imod=1;
991 Float_t factor = 1, frac = 0;
995 //Loop on the cells, get the cell amplitude and recalibration factor, multiply and and to the new energy
996 for(Int_t icell = 0; icell < ncells; icell++){
997 absId = index[icell];
998 frac = fraction[icell];
999 if(frac < 1e-5) frac = 1; //in case of EMCAL, this is set as 0 since unfolding is off
1001 if(!fCellsRecalibrated && IsRecalibrationOn()) {
1003 Int_t iTower = -1, iIphi = -1, iIeta = -1;
1004 geom->GetCellIndex(absId,imod,iTower,iIphi,iIeta);
1005 if(fEMCALRecalibrationFactors->GetEntries() <= imod) continue;
1006 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,irow,icol);
1007 factor = GetEMCALChannelRecalibrationFactor(imod,icol,irow);
1009 AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - recalibrate cell: module %d, col %d, row %d, cell fraction %f,recalibration factor %f, cell energy %f\n",
1010 imod,icol,irow,frac,factor,cells->GetCellAmplitude(absId)));
1014 energy += cells->GetCellAmplitude(absId)*factor*frac;
1016 if(emax < cells->GetCellAmplitude(absId)*factor*frac){
1017 emax = cells->GetCellAmplitude(absId)*factor*frac;
1022 cluster->SetE(energy);
1024 AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - Energy before %f, after %f\n",cluster->E(),energy));
1026 // Recalculate time of cluster only for ESDs
1027 if(!strcmp("AliESDCaloCluster",Form("%s",cluster->ClassName()))) {
1030 Double_t weightedTime = 0;
1031 Double_t weight = 0;
1032 Double_t weightTot = 0;
1033 Double_t maxcellTime = 0;
1034 for(Int_t icell = 0; icell < ncells; icell++){
1035 absId = index[icell];
1036 frac = fraction[icell];
1037 if(frac < 1e-5) frac = 1; //in case of EMCAL, this is set as 0 since unfolding is off
1039 Double_t celltime = cells->GetCellTime(absId);
1040 RecalibrateCellTime(absId, bc, celltime);
1041 if(absId == absIdMax) maxcellTime = celltime;
1043 if(!fCellsRecalibrated){
1045 Int_t iTower = -1, iIphi = -1, iIeta = -1;
1046 geom->GetCellIndex(absId,imod,iTower,iIphi,iIeta);
1047 if(fEMCALRecalibrationFactors->GetEntries() <= imod) continue;
1048 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,irow,icol);
1049 factor = GetEMCALChannelRecalibrationFactor(imod,icol,irow);
1051 AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - recalibrate cell:"
1052 " module %d, col %d, row %d, cell fraction %f,recalibration factor %f, cell energy %f\n",
1053 imod,icol,irow,frac,factor,cells->GetCellTime(absId)));
1057 weight = GetCellWeight(cells->GetCellAmplitude(absId)*factor*frac , energy );
1058 weightTot += weight;
1059 weightedTime += celltime * weight;
1063 cluster->SetTOF(weightedTime/weightTot);
1065 cluster->SetTOF(maxcellTime);
1069 //_____________________________________________________________
1070 void AliEMCALRecoUtils::RecalibrateCells(AliVCaloCells * cells,
1073 // Recalibrate the cells time and energy, considering the recalibration map and the energy
1074 // of the cells that compose the cluster.
1075 // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1077 if(!IsRecalibrationOn() && !IsTimeRecalibrationOn()) return;
1080 AliInfo("Cells pointer null!");
1085 Bool_t accept = kFALSE;
1089 Int_t nEMcell = cells->GetNumberOfCells() ;
1090 for (Int_t iCell = 0; iCell < nEMcell; iCell++) {
1092 absId = cells->GetCellNumber(iCell);
1094 accept = AcceptCalibrateCell(absId, bc, ecell ,tcell ,cells);
1101 cells->SetCell(iCell,absId,ecell, tcell);
1104 fCellsRecalibrated = kTRUE;
1107 //_______________________________________________________________________________________________________
1108 void AliEMCALRecoUtils::RecalibrateCellTime(const Int_t absId, const Int_t bc, Double_t & celltime) const
1110 // Recalibrate time of cell with absID considering the recalibration map
1111 // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1113 if(!fCellsRecalibrated && IsTimeRecalibrationOn() && bc >= 0){
1114 celltime -= GetEMCALChannelTimeRecalibrationFactor(bc%4,absId)*1.e-9; ;
1118 //______________________________________________________________________________
1119 void AliEMCALRecoUtils::RecalculateClusterPosition(const AliEMCALGeometry *geom,
1120 AliVCaloCells* cells,
1123 //For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1126 AliInfo("Cluster pointer null!");
1130 if (fPosAlgo==kPosTowerGlobal) RecalculateClusterPositionFromTowerGlobal( geom, cells, clu);
1131 else if(fPosAlgo==kPosTowerIndex) RecalculateClusterPositionFromTowerIndex ( geom, cells, clu);
1132 else AliDebug(2,"Algorithm to recalculate position not selected, do nothing.");
1135 //_____________________________________________________________________________________________
1136 void AliEMCALRecoUtils::RecalculateClusterPositionFromTowerGlobal(const AliEMCALGeometry *geom,
1137 AliVCaloCells* cells,
1140 // For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1141 // The algorithm is a copy of what is done in AliEMCALRecPoint
1143 Double_t eCell = 0.;
1144 Float_t fraction = 1.;
1145 Float_t recalFactor = 1.;
1148 Int_t iTower = -1, iIphi = -1, iIeta = -1;
1149 Int_t iSupModMax = -1, iSM=-1, iphi = -1, ieta = -1;
1150 Float_t weight = 0., totalWeight=0.;
1151 Float_t newPos[3] = {0,0,0};
1152 Double_t pLocal[3], pGlobal[3];
1153 Bool_t shared = kFALSE;
1155 Float_t clEnergy = clu->E(); //Energy already recalibrated previously
1156 GetMaxEnergyCell(geom, cells, clu, absId, iSupModMax, ieta, iphi,shared);
1157 Double_t depth = GetDepth(clEnergy,fParticleType,iSupModMax) ;
1159 //printf("** Cluster energy %f, ncells %d, depth %f\n",clEnergy,clu->GetNCells(),depth);
1161 for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++) {
1163 absId = clu->GetCellAbsId(iDig);
1164 fraction = clu->GetCellAmplitudeFraction(iDig);
1165 if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1167 if (!fCellsRecalibrated){
1168 geom->GetCellIndex(absId,iSM,iTower,iIphi,iIeta);
1169 geom->GetCellPhiEtaIndexInSModule(iSM,iTower,iIphi, iIeta,iphi,ieta);
1171 if(IsRecalibrationOn()) {
1172 recalFactor = GetEMCALChannelRecalibrationFactor(iSM,ieta,iphi);
1176 eCell = cells->GetCellAmplitude(absId)*fraction*recalFactor;
1178 weight = GetCellWeight(eCell,clEnergy);
1179 totalWeight += weight;
1181 geom->RelPosCellInSModule(absId,depth,pLocal[0],pLocal[1],pLocal[2]);
1182 //printf("pLocal (%f,%f,%f), SM %d, absId %d\n",pLocal[0],pLocal[1],pLocal[2],iSupModMax,absId);
1183 geom->GetGlobal(pLocal,pGlobal,iSupModMax);
1184 //printf("pLocal (%f,%f,%f)\n",pGlobal[0],pGlobal[1],pGlobal[2]);
1186 for(int i=0; i<3; i++ ) newPos[i] += (weight*pGlobal[i]);
1190 for(int i=0; i<3; i++ ) newPos[i] /= totalWeight;
1193 //Float_t pos[]={0,0,0};
1194 //clu->GetPosition(pos);
1195 //printf("OldPos : %2.3f,%2.3f,%2.3f\n",pos[0],pos[1],pos[2]);
1196 //printf("NewPos : %2.3f,%2.3f,%2.3f\n",newPos[0],newPos[1],newPos[2]);
1198 if(iSupModMax > 1) {//sector 1
1199 newPos[0] +=fMisalTransShift[3];//-=3.093;
1200 newPos[1] +=fMisalTransShift[4];//+=6.82;
1201 newPos[2] +=fMisalTransShift[5];//+=1.635;
1202 //printf(" + : %2.3f,%2.3f,%2.3f\n",fMisalTransShift[3],fMisalTransShift[4],fMisalTransShift[5]);
1204 newPos[0] +=fMisalTransShift[0];//+=1.134;
1205 newPos[1] +=fMisalTransShift[1];//+=8.2;
1206 newPos[2] +=fMisalTransShift[2];//+=1.197;
1207 //printf(" + : %2.3f,%2.3f,%2.3f\n",fMisalTransShift[0],fMisalTransShift[1],fMisalTransShift[2]);
1209 //printf("NewPos : %2.3f,%2.3f,%2.3f\n",newPos[0],newPos[1],newPos[2]);
1211 clu->SetPosition(newPos);
1214 //____________________________________________________________________________________________
1215 void AliEMCALRecoUtils::RecalculateClusterPositionFromTowerIndex(const AliEMCALGeometry *geom,
1216 AliVCaloCells* cells,
1219 // For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1220 // The algorithm works with the tower indeces, averages the indeces and from them it calculates the global position
1222 Double_t eCell = 1.;
1223 Float_t fraction = 1.;
1224 Float_t recalFactor = 1.;
1228 Int_t iIphi = -1, iIeta = -1;
1229 Int_t iSupMod = -1, iSupModMax = -1;
1230 Int_t iphi = -1, ieta =-1;
1231 Bool_t shared = kFALSE;
1233 Float_t clEnergy = clu->E(); //Energy already recalibrated previously.
1234 GetMaxEnergyCell(geom, cells, clu, absId, iSupModMax, ieta, iphi,shared);
1235 Float_t depth = GetDepth(clEnergy,fParticleType,iSupMod) ;
1237 Float_t weight = 0., weightedCol = 0., weightedRow = 0., totalWeight=0.;
1238 Bool_t areInSameSM = kTRUE; //exclude clusters with cells in different SMs for now
1239 Int_t startingSM = -1;
1241 for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++) {
1242 absId = clu->GetCellAbsId(iDig);
1243 fraction = clu->GetCellAmplitudeFraction(iDig);
1244 if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1246 if (iDig==0) startingSM = iSupMod;
1247 else if(iSupMod != startingSM) areInSameSM = kFALSE;
1249 eCell = cells->GetCellAmplitude(absId);
1251 geom->GetCellIndex(absId,iSupMod,iTower,iIphi,iIeta);
1252 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);
1254 if (!fCellsRecalibrated){
1255 if(IsRecalibrationOn()) {
1256 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1260 eCell = cells->GetCellAmplitude(absId)*fraction*recalFactor;
1262 weight = GetCellWeight(eCell,clEnergy);
1263 if(weight < 0) weight = 0;
1264 totalWeight += weight;
1265 weightedCol += ieta*weight;
1266 weightedRow += iphi*weight;
1268 //printf("Max cell? cell %d, amplitude org %f, fraction %f, recalibration %f, amplitude new %f \n",cellAbsId, cells->GetCellAmplitude(cellAbsId), fraction, recalFactor, eCell) ;
1271 Float_t xyzNew[]={0.,0.,0.};
1272 if(areInSameSM == kTRUE) {
1273 //printf("In Same SM\n");
1274 weightedCol = weightedCol/totalWeight;
1275 weightedRow = weightedRow/totalWeight;
1276 geom->RecalculateTowerPosition(weightedRow, weightedCol, iSupModMax, depth, fMisalTransShift, fMisalRotShift, xyzNew);
1278 //printf("In Different SM\n");
1279 geom->RecalculateTowerPosition(iphi, ieta, iSupModMax, depth, fMisalTransShift, fMisalRotShift, xyzNew);
1282 clu->SetPosition(xyzNew);
1285 //___________________________________________________________________________________________
1286 void AliEMCALRecoUtils::RecalculateClusterDistanceToBadChannel(const AliEMCALGeometry * geom,
1287 AliVCaloCells* cells,
1288 AliVCluster * cluster)
1290 //re-evaluate distance to bad channel with updated bad map
1292 if(!fRecalDistToBadChannels) return;
1295 AliInfo("Cluster pointer null!");
1299 //Get channels map of the supermodule where the cluster is.
1300 Int_t absIdMax = -1, iSupMod =-1, icolM = -1, irowM = -1;
1301 Bool_t shared = kFALSE;
1302 GetMaxEnergyCell(geom, cells, cluster, absIdMax, iSupMod, icolM, irowM, shared);
1303 TH2D* hMap = (TH2D*)fEMCALBadChannelMap->At(iSupMod);
1306 Float_t minDist = 10000.;
1309 //Loop on tower status map
1310 for(Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++){
1311 for(Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++){
1312 //Check if tower is bad.
1313 if(hMap->GetBinContent(icol,irow)==0) continue;
1314 //printf("AliEMCALRecoUtils::RecalculateDistanceToBadChannels() - \n \t Bad channel in SM %d, col %d, row %d, \n \t Cluster max in col %d, row %d\n",
1315 // iSupMod,icol, irow, icolM,irowM);
1317 dRrow=TMath::Abs(irowM-irow);
1318 dRcol=TMath::Abs(icolM-icol);
1319 dist=TMath::Sqrt(dRrow*dRrow+dRcol*dRcol);
1321 //printf("MIN DISTANCE TO BAD %2.2f\n",dist);
1327 //In case the cluster is shared by 2 SuperModules, need to check the map of the second Super Module
1330 Int_t iSupMod2 = -1;
1332 //The only possible combinations are (0,1), (2,3) ... (8,9)
1333 if(iSupMod%2) iSupMod2 = iSupMod-1;
1334 else iSupMod2 = iSupMod+1;
1335 hMap2 = (TH2D*)fEMCALBadChannelMap->At(iSupMod2);
1337 //Loop on tower status map of second super module
1338 for(Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++){
1339 for(Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++){
1340 //Check if tower is bad.
1341 if(hMap2->GetBinContent(icol,irow)==0) continue;
1342 //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",
1343 // iSupMod2,icol, irow,iSupMod,icolM,irowM);
1344 dRrow=TMath::Abs(irow-irowM);
1347 dRcol=TMath::Abs(icol-(AliEMCALGeoParams::fgkEMCALCols+icolM));
1349 dRcol=TMath::Abs(AliEMCALGeoParams::fgkEMCALCols+icol-icolM);
1352 dist=TMath::Sqrt(dRrow*dRrow+dRcol*dRcol);
1353 if(dist < minDist) minDist = dist;
1356 }// shared cluster in 2 SuperModules
1358 AliDebug(2,Form("Max cluster cell (SM,col,row)=(%d %d %d) - Distance to Bad Channel %2.2f",iSupMod, icolM, irowM, minDist));
1359 cluster->SetDistanceToBadChannel(minDist);
1362 //__________________________________________________________________
1363 void AliEMCALRecoUtils::RecalculateClusterPID(AliVCluster * cluster)
1365 //re-evaluate identification parameters with bayesian
1368 AliInfo("Cluster pointer null!");
1372 if ( cluster->GetM02() != 0)
1373 fPIDUtils->ComputePID(cluster->E(),cluster->GetM02());
1375 Float_t pidlist[AliPID::kSPECIESN+1];
1376 for(Int_t i = 0; i < AliPID::kSPECIESN+1; i++) pidlist[i] = fPIDUtils->GetPIDFinal(i);
1378 cluster->SetPID(pidlist);
1381 //____________________________________________________________________________________________
1382 void AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(const AliEMCALGeometry * geom,
1383 AliVCaloCells* cells,
1384 AliVCluster * cluster)
1386 // Calculates new center of gravity in the local EMCAL-module coordinates
1387 // and tranfers into global ALICE coordinates
1388 // Calculates Dispersion and main axis
1391 AliInfo("Cluster pointer null!");
1397 Double_t eCell = 0.;
1398 Float_t fraction = 1.;
1399 Float_t recalFactor = 1.;
1407 Double_t etai = -1.;
1408 Double_t phii = -1.;
1415 Double_t xmean = 0.;
1416 Double_t zmean = 0.;
1419 for(Int_t iDigit=0; iDigit < cluster->GetNCells(); iDigit++) {
1421 //Get from the absid the supermodule, tower and eta/phi numbers
1422 geom->GetCellIndex(cluster->GetCellAbsId(iDigit),iSupMod,iTower,iIphi,iIeta);
1423 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi,iIeta, iphi,ieta);
1425 //Get the cell energy, if recalibration is on, apply factors
1426 fraction = cluster->GetCellAmplitudeFraction(iDigit);
1427 if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1429 if (!fCellsRecalibrated){
1430 if(IsRecalibrationOn()) {
1431 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1435 eCell = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1437 if(cluster->E() > 0 && eCell > 0){
1439 w = GetCellWeight(eCell,cluster->E());
1441 etai=(Double_t)ieta;
1442 phii=(Double_t)iphi;
1447 dxx += w * etai * etai ;
1449 dzz += w * phii * phii ;
1451 dxz += w * etai * phii ;
1454 AliError(Form("Wrong energy %f and/or amplitude %f\n", eCell, cluster->E()));
1457 //Normalize to the weight
1463 AliError(Form("Wrong weight %f\n", wtot));
1465 //Calculate dispersion
1466 for(Int_t iDigit=0; iDigit < cluster->GetNCells(); iDigit++) {
1468 //Get from the absid the supermodule, tower and eta/phi numbers
1469 geom->GetCellIndex(cluster->GetCellAbsId(iDigit),iSupMod,iTower,iIphi,iIeta);
1470 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi,iIeta, iphi,ieta);
1472 //Get the cell energy, if recalibration is on, apply factors
1473 fraction = cluster->GetCellAmplitudeFraction(iDigit);
1474 if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1475 if (IsRecalibrationOn()) {
1476 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1478 eCell = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1480 if(cluster->E() > 0 && eCell > 0){
1482 w = GetCellWeight(eCell,cluster->E());
1484 etai=(Double_t)ieta;
1485 phii=(Double_t)iphi;
1486 if(w > 0.0) d += w*((etai-xmean)*(etai-xmean)+(phii-zmean)*(phii-zmean));
1489 AliError(Form("Wrong energy %f and/or amplitude %f\n", eCell, cluster->E()));
1492 //Normalize to the weigth and set shower shape parameters
1493 if (wtot > 0 && nstat > 1) {
1498 dxx -= xmean * xmean ;
1499 dzz -= zmean * zmean ;
1500 dxz -= xmean * zmean ;
1501 cluster->SetM02(0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ));
1502 cluster->SetM20(0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ));
1506 cluster->SetM20(0.) ;
1507 cluster->SetM02(0.) ;
1511 cluster->SetDispersion(TMath::Sqrt(d)) ;
1513 cluster->SetDispersion(0) ;
1516 //____________________________________________________________________________
1517 void AliEMCALRecoUtils::FindMatches(AliVEvent *event,
1518 TObjArray * clusterArr,
1519 const AliEMCALGeometry *geom)
1521 //This function should be called before the cluster loop
1522 //Before call this function, please recalculate the cluster positions
1523 //Given the input event, loop over all the tracks, select the closest cluster as matched with fCutR
1524 //Store matched cluster indexes and residuals
1526 fMatchedTrackIndex ->Reset();
1527 fMatchedClusterIndex->Reset();
1528 fResidualPhi->Reset();
1529 fResidualEta->Reset();
1531 fMatchedTrackIndex ->Set(1000);
1532 fMatchedClusterIndex->Set(1000);
1533 fResidualPhi->Set(1000);
1534 fResidualEta->Set(1000);
1536 AliESDEvent* esdevent = dynamic_cast<AliESDEvent*> (event);
1537 AliAODEvent* aodevent = dynamic_cast<AliAODEvent*> (event);
1539 // init the magnetic field if not already on
1540 if(!TGeoGlobalMagField::Instance()->GetField()){
1541 AliInfo("Init the magnetic field\n");
1544 esdevent->InitMagneticField();
1548 Double_t curSol = 30000*aodevent->GetMagneticField()/5.00668;
1549 Double_t curDip = 6000 *aodevent->GetMuonMagFieldScale();
1550 AliMagF *field = AliMagF::CreateFieldMap(curSol,curDip);
1551 TGeoGlobalMagField::Instance()->SetField(field);
1555 AliInfo("Mag Field not initialized, null esd/aod evetn pointers");
1560 TObjArray *clusterArray = 0x0;
1563 clusterArray = new TObjArray(event->GetNumberOfCaloClusters());
1564 for(Int_t icl=0; icl<event->GetNumberOfCaloClusters(); icl++)
1566 AliVCluster *cluster = (AliVCluster*) event->GetCaloCluster(icl);
1567 if(geom && !IsGoodCluster(cluster,geom,(AliVCaloCells*)event->GetEMCALCells())) continue;
1568 clusterArray->AddAt(cluster,icl);
1574 for (Int_t i=0; i<21;i++) cv[i]=0;
1575 for(Int_t itr=0; itr<event->GetNumberOfTracks(); itr++)
1577 AliExternalTrackParam *trackParam = 0;
1579 //If the input event is ESD, the starting point for extrapolation is TPCOut, if available, or TPCInner
1580 AliESDtrack *esdTrack = 0;
1581 AliAODTrack *aodTrack = 0;
1584 esdTrack = esdevent->GetTrack(itr);
1585 if(!esdTrack) continue;
1586 if(!IsAccepted(esdTrack)) continue;
1587 if(esdTrack->Pt()<fCutMinTrackPt) continue;
1588 Double_t phi = esdTrack->Phi()*TMath::RadToDeg();
1589 if(TMath::Abs(esdTrack->Eta())>0.8 || phi <= 20 || phi >= 240 ) continue;
1590 trackParam = const_cast<AliExternalTrackParam*>(esdTrack->GetInnerParam());
1593 //If the input event is AOD, the starting point for extrapolation is at vertex
1594 //AOD tracks are selected according to its filterbit.
1597 aodTrack = aodevent->GetTrack(itr);
1598 if(!aodTrack) continue;
1599 if(!aodTrack->TestFilterMask(fAODFilterMask)) continue; //Select AOD tracks that fulfill GetStandardITSTPCTrackCuts2010()
1600 if(aodTrack->Pt()<fCutMinTrackPt) continue;
1601 Double_t phi = aodTrack->Phi()*TMath::RadToDeg();
1602 if(TMath::Abs(aodTrack->Eta())>0.8 || phi <= 20 || phi >= 240 ) continue;
1603 Double_t pos[3],mom[3];
1604 aodTrack->GetXYZ(pos);
1605 aodTrack->GetPxPyPz(mom);
1606 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()));
1607 trackParam= new AliExternalTrackParam(pos,mom,cv,aodTrack->Charge());
1610 //Return if the input data is not "AOD" or "ESD"
1613 printf("Wrong input data type! Should be \"AOD\" or \"ESD\"\n");
1616 clusterArray->Clear();
1617 delete clusterArray;
1622 if(!trackParam) continue;
1624 //Extrapolate the track to EMCal surface
1625 AliExternalTrackParam emcalParam(*trackParam);
1627 if(!ExtrapolateTrackToEMCalSurface(&emcalParam, 430., fMass, fStepSurface, eta, phi))
1629 if(aodevent && trackParam) delete trackParam;
1635 // esdTrack->SetOuterParam(&emcalParam,AliExternalTrackParam::kMultSec);
1638 if(TMath::Abs(eta)>0.75 || (phi) < 70*TMath::DegToRad() || (phi) > 190*TMath::DegToRad())
1640 if(aodevent && trackParam) delete trackParam;
1645 //Find matched clusters
1647 Float_t dEta = -999, dPhi = -999;
1650 index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArray, dEta, dPhi);
1654 index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArr, dEta, dPhi);
1659 fMatchedTrackIndex ->AddAt(itr,matched);
1660 fMatchedClusterIndex ->AddAt(index,matched);
1661 fResidualEta ->AddAt(dEta,matched);
1662 fResidualPhi ->AddAt(dPhi,matched);
1665 if(aodevent && trackParam) delete trackParam;
1670 clusterArray->Clear();
1671 delete clusterArray;
1674 AliDebug(2,Form("Number of matched pairs = %d !\n",matched));
1676 fMatchedTrackIndex ->Set(matched);
1677 fMatchedClusterIndex ->Set(matched);
1678 fResidualPhi ->Set(matched);
1679 fResidualEta ->Set(matched);
1682 //________________________________________________________________________________
1683 Int_t AliEMCALRecoUtils::FindMatchedClusterInEvent(const AliESDtrack *track,
1684 const AliVEvent *event,
1685 const AliEMCALGeometry *geom,
1686 Float_t &dEta, Float_t &dPhi)
1689 // This function returns the index of matched cluster to input track
1690 // Returns -1 if no match is found
1692 Double_t phiV = track->Phi()*TMath::RadToDeg();
1693 if(TMath::Abs(track->Eta())>0.8 || phiV <= 20 || phiV >= 240 ) return index;
1694 AliExternalTrackParam *trackParam = const_cast<AliExternalTrackParam*>(track->GetInnerParam());
1695 if(!trackParam) return index;
1696 AliExternalTrackParam emcalParam(*trackParam);
1698 if(!ExtrapolateTrackToEMCalSurface(&emcalParam, 430., fMass, fStepSurface, eta, phi)) return index;
1699 if(TMath::Abs(eta)>0.75 || (phi) < 70*TMath::DegToRad() || (phi) > 190*TMath::DegToRad()) return index;
1701 TObjArray *clusterArr = new TObjArray(event->GetNumberOfCaloClusters());
1703 for(Int_t icl=0; icl<event->GetNumberOfCaloClusters(); icl++)
1705 AliVCluster *cluster = (AliVCluster*) event->GetCaloCluster(icl);
1706 if(geom && !IsGoodCluster(cluster,geom,(AliVCaloCells*)event->GetEMCALCells())) continue;
1707 clusterArr->AddAt(cluster,icl);
1710 index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArr, dEta, dPhi);
1711 clusterArr->Clear();
1717 //________________________________________________________________________________________
1718 Int_t AliEMCALRecoUtils::FindMatchedClusterInClusterArr(AliExternalTrackParam *emcalParam,
1719 AliExternalTrackParam *trkParam,
1720 TObjArray * clusterArr,
1721 Float_t &dEta, Float_t &dPhi)
1723 dEta=-999, dPhi=-999;
1724 Float_t dRMax = fCutR, dEtaMax=fCutEta, dPhiMax=fCutPhi;
1726 Float_t tmpEta=-999, tmpPhi=-999;
1728 Double_t exPos[3] = {0.,0.,0.};
1729 if(!emcalParam->GetXYZ(exPos)) return index;
1731 Float_t clsPos[3] = {0.,0.,0.};
1732 for(Int_t icl=0; icl<clusterArr->GetEntriesFast(); icl++)
1734 AliVCluster *cluster = dynamic_cast<AliVCluster*> (clusterArr->At(icl)) ;
1735 if(!cluster || !cluster->IsEMCAL()) continue;
1736 cluster->GetPosition(clsPos);
1737 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));
1738 if(dR > fClusterWindow) continue;
1740 AliExternalTrackParam trkPamTmp (*trkParam);//Retrieve the starting point every time before the extrapolation
1741 if(!ExtrapolateTrackToCluster(&trkPamTmp, cluster, fMass, fStepCluster, tmpEta, tmpPhi)) continue;
1744 Float_t tmpR=TMath::Sqrt(tmpEta*tmpEta + tmpPhi*tmpPhi);
1753 else if(fCutEtaPhiSeparate)
1755 if(TMath::Abs(tmpEta)<TMath::Abs(dEtaMax) && TMath::Abs(tmpPhi)<TMath::Abs(dPhiMax))
1764 printf("Error: please specify your cut criteria\n");
1765 printf("To cut on sqrt(dEta^2+dPhi^2), use: SwitchOnCutEtaPhiSum()\n");
1766 printf("To cut on dEta and dPhi separately, use: SwitchOnCutEtaPhiSeparate()\n");
1777 //------------------------------------------------------------------------------------
1778 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToEMCalSurface(AliExternalTrackParam *trkParam,
1779 const Double_t emcalR,
1780 const Double_t mass,
1781 const Double_t step,
1785 //Extrapolate track to EMCAL surface
1787 eta = -999, phi = -999;
1788 if(!trkParam) return kFALSE;
1789 if(!AliTrackerBase::PropagateTrackToBxByBz(trkParam, emcalR, mass, step, kTRUE, 0.8, -1)) return kFALSE;
1790 Double_t trkPos[3] = {0.,0.,0.};
1791 if(!trkParam->GetXYZ(trkPos)) return kFALSE;
1792 TVector3 trkPosVec(trkPos[0],trkPos[1],trkPos[2]);
1793 eta = trkPosVec.Eta();
1794 phi = trkPosVec.Phi();
1796 phi += 2*TMath::Pi();
1801 //-----------------------------------------------------------------------------------
1802 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToPosition(AliExternalTrackParam *trkParam,
1803 const Float_t *clsPos,
1810 //Return the residual by extrapolating a track param to a global position
1814 if(!trkParam) return kFALSE;
1815 Double_t trkPos[3] = {0.,0.,0.};
1816 TVector3 vec(clsPos[0],clsPos[1],clsPos[2]);
1817 Double_t alpha = ((int)(vec.Phi()*TMath::RadToDeg()/20)+0.5)*20*TMath::DegToRad();
1818 vec.RotateZ(-alpha); //Rotate the cluster to the local extrapolation coordinate system
1819 if(!AliTrackerBase::PropagateTrackToBxByBz(trkParam, vec.X(), mass, step,kTRUE, 0.8, -1)) return kFALSE;
1820 if(!trkParam->GetXYZ(trkPos)) return kFALSE; //Get the extrapolated global position
1822 TVector3 clsPosVec(clsPos[0],clsPos[1],clsPos[2]);
1823 TVector3 trkPosVec(trkPos[0],trkPos[1],trkPos[2]);
1825 // track cluster matching
1826 tmpPhi = clsPosVec.DeltaPhi(trkPosVec); // tmpPhi is between -pi and pi
1827 tmpEta = clsPosVec.Eta()-trkPosVec.Eta();
1832 //----------------------------------------------------------------------------------
1833 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToCluster(AliExternalTrackParam *trkParam,
1834 AliVCluster *cluster,
1835 const Double_t mass,
1836 const Double_t step,
1841 //Return the residual by extrapolating a track param to a cluster
1845 if(!cluster || !trkParam) return kFALSE;
1847 Float_t clsPos[3] = {0.,0.,0.};
1848 cluster->GetPosition(clsPos);
1850 return ExtrapolateTrackToPosition(trkParam, clsPos, mass, step, tmpEta, tmpPhi);
1853 //---------------------------------------------------------------------------------
1854 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToCluster(AliExternalTrackParam *trkParam,
1855 AliVCluster *cluster,
1860 //Return the residual by extrapolating a track param to a clusterfStepCluster
1863 return ExtrapolateTrackToCluster(trkParam, cluster, fMass, fStepCluster, tmpEta, tmpPhi);
1866 //_______________________________________________________________________
1867 void AliEMCALRecoUtils::GetMatchedResiduals(const Int_t clsIndex,
1868 Float_t &dEta, Float_t &dPhi)
1870 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
1871 //Get the residuals dEta and dPhi for this cluster to the closest track
1872 //Works with ESDs and AODs
1874 if( FindMatchedPosForCluster(clsIndex) >= 999 )
1876 AliDebug(2,"No matched tracks found!\n");
1881 dEta = fResidualEta->At(FindMatchedPosForCluster(clsIndex));
1882 dPhi = fResidualPhi->At(FindMatchedPosForCluster(clsIndex));
1885 //______________________________________________________________________________________________
1886 void AliEMCALRecoUtils::GetMatchedClusterResiduals(Int_t trkIndex, Float_t &dEta, Float_t &dPhi)
1888 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
1889 //Get the residuals dEta and dPhi for this track to the closest cluster
1890 //Works with ESDs and AODs
1892 if( FindMatchedPosForTrack(trkIndex) >= 999 )
1894 AliDebug(2,"No matched cluster found!\n");
1899 dEta = fResidualEta->At(FindMatchedPosForTrack(trkIndex));
1900 dPhi = fResidualPhi->At(FindMatchedPosForTrack(trkIndex));
1903 //__________________________________________________________
1904 Int_t AliEMCALRecoUtils::GetMatchedTrackIndex(Int_t clsIndex)
1906 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
1907 //Get the index of matched track to this cluster
1908 //Works with ESDs and AODs
1910 if(IsClusterMatched(clsIndex))
1911 return fMatchedTrackIndex->At(FindMatchedPosForCluster(clsIndex));
1916 //__________________________________________________________
1917 Int_t AliEMCALRecoUtils::GetMatchedClusterIndex(Int_t trkIndex)
1919 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
1920 //Get the index of matched cluster to this track
1921 //Works with ESDs and AODs
1923 if(IsTrackMatched(trkIndex))
1924 return fMatchedClusterIndex->At(FindMatchedPosForTrack(trkIndex));
1929 //__________________________________________________
1930 Bool_t AliEMCALRecoUtils::IsClusterMatched(Int_t clsIndex) const
1932 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
1933 //Returns if the cluster has a match
1934 if(FindMatchedPosForCluster(clsIndex) < 999)
1940 //__________________________________________________
1941 Bool_t AliEMCALRecoUtils::IsTrackMatched(Int_t trkIndex) const
1943 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
1944 //Returns if the track has a match
1945 if(FindMatchedPosForTrack(trkIndex) < 999)
1951 //__________________________________________________________
1952 UInt_t AliEMCALRecoUtils::FindMatchedPosForCluster(Int_t clsIndex) const
1954 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
1955 //Returns the position of the match in the fMatchedClusterIndex array
1956 Float_t tmpR = fCutR;
1959 for(Int_t i=0; i<fMatchedClusterIndex->GetSize(); i++) {
1960 if(fMatchedClusterIndex->At(i)==clsIndex) {
1961 Float_t r = TMath::Sqrt(fResidualEta->At(i)*fResidualEta->At(i) + fResidualPhi->At(i)*fResidualPhi->At(i));
1965 AliDebug(3,Form("Matched cluster index: index: %d, dEta: %2.4f, dPhi: %2.4f.\n",fMatchedClusterIndex->At(i),fResidualEta->At(i),fResidualPhi->At(i)));
1972 //__________________________________________________________
1973 UInt_t AliEMCALRecoUtils::FindMatchedPosForTrack(Int_t trkIndex) const
1975 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
1976 //Returns the position of the match in the fMatchedTrackIndex array
1977 Float_t tmpR = fCutR;
1980 for(Int_t i=0; i<fMatchedTrackIndex->GetSize(); i++) {
1981 if(fMatchedTrackIndex->At(i)==trkIndex) {
1982 Float_t r = TMath::Sqrt(fResidualEta->At(i)*fResidualEta->At(i) + fResidualPhi->At(i)*fResidualPhi->At(i));
1986 AliDebug(3,Form("Matched track index: index: %d, dEta: %2.4f, dPhi: %2.4f.\n",fMatchedTrackIndex->At(i),fResidualEta->At(i),fResidualPhi->At(i)));
1993 //__________________________________________________________________________
1994 Bool_t AliEMCALRecoUtils::IsGoodCluster(AliVCluster *cluster,
1995 const AliEMCALGeometry *geom,
1996 AliVCaloCells* cells,const Int_t bc)
1998 // check if the cluster survives some quality cut
2001 Bool_t isGood=kTRUE;
2003 if(!cluster || !cluster->IsEMCAL()) return kFALSE;
2005 if(ClusterContainsBadChannel(geom,cluster->GetCellsAbsId(),cluster->GetNCells())) return kFALSE;
2007 if(!CheckCellFiducialRegion(geom,cluster,cells)) return kFALSE;
2009 if(IsExoticCluster(cluster, cells,bc)) return kFALSE;
2014 //__________________________________________________________
2015 Bool_t AliEMCALRecoUtils::IsAccepted(AliESDtrack *esdTrack)
2017 // Given a esd track, return whether the track survive all the cuts
2019 // The different quality parameter are first
2020 // retrieved from the track. then it is found out what cuts the
2021 // track did not survive and finally the cuts are imposed.
2023 UInt_t status = esdTrack->GetStatus();
2025 Int_t nClustersITS = esdTrack->GetITSclusters(0);
2026 Int_t nClustersTPC = esdTrack->GetTPCclusters(0);
2028 Float_t chi2PerClusterITS = -1;
2029 Float_t chi2PerClusterTPC = -1;
2030 if (nClustersITS!=0)
2031 chi2PerClusterITS = esdTrack->GetITSchi2()/Float_t(nClustersITS);
2032 if (nClustersTPC!=0)
2033 chi2PerClusterTPC = esdTrack->GetTPCchi2()/Float_t(nClustersTPC);
2037 if(fTrackCutsType==kGlobalCut)
2039 Float_t maxDCAToVertexXYPtDep = 0.0182 + 0.0350/TMath::Power(esdTrack->Pt(),1.01); //This expression comes from AliESDtrackCuts::GetStandardITSTPCTrackCuts2010()
2040 //AliDebug(3,Form("Track pT = %f, DCAtoVertexXY = %f",esdTrack->Pt(),MaxDCAToVertexXYPtDep));
2041 SetMaxDCAToVertexXY(maxDCAToVertexXYPtDep); //Set pT dependent DCA cut to vertex in x-y plane
2047 esdTrack->GetImpactParameters(b,bCov);
2048 if (bCov[0]<=0 || bCov[2]<=0) {
2049 AliDebug(1, "Estimated b resolution lower or equal zero!");
2050 bCov[0]=0; bCov[2]=0;
2053 Float_t dcaToVertexXY = b[0];
2054 Float_t dcaToVertexZ = b[1];
2055 Float_t dcaToVertex = -1;
2057 if (fCutDCAToVertex2D)
2058 dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY/fCutMaxDCAToVertexXY/fCutMaxDCAToVertexXY + dcaToVertexZ*dcaToVertexZ/fCutMaxDCAToVertexZ/fCutMaxDCAToVertexZ);
2060 dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY + dcaToVertexZ*dcaToVertexZ);
2064 Bool_t cuts[kNCuts];
2065 for (Int_t i=0; i<kNCuts; i++) cuts[i]=kFALSE;
2067 // track quality cuts
2068 if (fCutRequireTPCRefit && (status&AliESDtrack::kTPCrefit)==0)
2070 if (fCutRequireITSRefit && (status&AliESDtrack::kITSrefit)==0)
2072 if (nClustersTPC<fCutMinNClusterTPC)
2074 if (nClustersITS<fCutMinNClusterITS)
2076 if (chi2PerClusterTPC>fCutMaxChi2PerClusterTPC)
2078 if (chi2PerClusterITS>fCutMaxChi2PerClusterITS)
2080 if (!fCutAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0)
2082 if (fCutDCAToVertex2D && dcaToVertex > 1)
2084 if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexXY) > fCutMaxDCAToVertexXY)
2086 if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexZ) > fCutMaxDCAToVertexZ)
2089 if(fTrackCutsType==kGlobalCut)
2091 //Require at least one SPD point + anything else in ITS
2092 if( (esdTrack->HasPointOnITSLayer(0) || esdTrack->HasPointOnITSLayer(1)) == kFALSE)
2097 for (Int_t i=0; i<kNCuts; i++)
2098 if (cuts[i]) {cut = kTRUE;}
2107 //_____________________________________
2108 void AliEMCALRecoUtils::InitTrackCuts()
2110 //Intilize the track cut criteria
2111 //By default these cuts are set according to AliESDtrackCuts::GetStandardTPCOnlyTrackCuts()
2112 //Also you can customize the cuts using the setters
2114 switch (fTrackCutsType)
2118 AliInfo(Form("Track cuts for matching: GetStandardTPCOnlyTrackCuts()"));
2120 SetMinNClustersTPC(70);
2121 SetMaxChi2PerClusterTPC(4);
2122 SetAcceptKinkDaughters(kFALSE);
2123 SetRequireTPCRefit(kFALSE);
2126 SetRequireITSRefit(kFALSE);
2127 SetMaxDCAToVertexZ(3.2);
2128 SetMaxDCAToVertexXY(2.4);
2129 SetDCAToVertex2D(kTRUE);
2136 AliInfo(Form("Track cuts for matching: GetStandardITSTPCTrackCuts2010(kTURE)"));
2138 SetMinNClustersTPC(70);
2139 SetMaxChi2PerClusterTPC(4);
2140 SetAcceptKinkDaughters(kFALSE);
2141 SetRequireTPCRefit(kTRUE);
2144 SetRequireITSRefit(kTRUE);
2145 SetMaxDCAToVertexZ(2);
2146 SetMaxDCAToVertexXY();
2147 SetDCAToVertex2D(kFALSE);
2154 AliInfo(Form("Track cuts for matching: Loose cut w/o DCA cut"));
2155 SetMinNClustersTPC(50);
2156 SetAcceptKinkDaughters(kTRUE);
2164 //________________________________________________________________________
2165 void AliEMCALRecoUtils::SetClusterMatchedToTrack(const AliESDEvent *event)
2167 // Checks if tracks are matched to EMC clusters and set the matched EMCAL cluster index to ESD track.
2169 Int_t nTracks = event->GetNumberOfTracks();
2170 for (Int_t iTrack = 0; iTrack < nTracks; ++iTrack) {
2171 AliESDtrack* track = event->GetTrack(iTrack);
2173 AliWarning(Form("Could not receive track %d", iTrack));
2176 Int_t matchClusIndex = GetMatchedClusterIndex(iTrack);
2177 track->SetEMCALcluster(matchClusIndex); //sets -1 if track not matched within residual
2178 if(matchClusIndex != -1)
2179 track->SetStatus(AliESDtrack::kEMCALmatch);
2181 track->ResetStatus(AliESDtrack::kEMCALmatch);
2183 AliDebug(2,"Track matched to closest cluster");
2186 //_________________________________________________________________________
2187 void AliEMCALRecoUtils::SetTracksMatchedToCluster(const AliESDEvent *event)
2189 // Checks if EMC clusters are matched to ESD track.
2190 // Adds track indexes of all the tracks matched to a cluster withing residuals in ESDCalocluster.
2192 for (Int_t iClus=0; iClus < event->GetNumberOfCaloClusters(); ++iClus) {
2193 AliESDCaloCluster *cluster = event->GetCaloCluster(iClus);
2194 if (!cluster->IsEMCAL())
2197 Int_t nTracks = event->GetNumberOfTracks();
2198 TArrayI arrayTrackMatched(nTracks);
2200 // Get the closest track matched to the cluster
2202 Int_t matchTrackIndex = GetMatchedTrackIndex(iClus);
2203 if (matchTrackIndex != -1) {
2204 arrayTrackMatched[nMatched] = matchTrackIndex;
2208 // Get all other tracks matched to the cluster
2209 for(Int_t iTrk=0; iTrk<nTracks; ++iTrk) {
2210 AliESDtrack* track = event->GetTrack(iTrk);
2211 if(iTrk == matchTrackIndex) continue;
2212 if(track->GetEMCALcluster() == iClus){
2213 arrayTrackMatched[nMatched] = iTrk;
2218 //printf("Tender::SetTracksMatchedToCluster - cluster E %f, N matches %d, first match %d\n",cluster->E(),nMatched,arrayTrackMatched[0]);
2220 arrayTrackMatched.Set(nMatched);
2221 cluster->AddTracksMatched(arrayTrackMatched);
2223 Float_t eta= -999, phi = -999;
2224 if (matchTrackIndex != -1)
2225 GetMatchedResiduals(iClus, eta, phi);
2226 cluster->SetTrackDistance(phi, eta);
2229 AliDebug(2,"Cluster matched to tracks");
2233 //___________________________________________________
2234 void AliEMCALRecoUtils::Print(const Option_t *) const
2238 printf("AliEMCALRecoUtils Settings: \n");
2239 printf("Misalignment shifts\n");
2240 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,
2241 fMisalTransShift[i*3],fMisalTransShift[i*3+1],fMisalTransShift[i*3+2],
2242 fMisalRotShift[i*3], fMisalRotShift[i*3+1], fMisalRotShift[i*3+2] );
2243 printf("Non linearity function %d, parameters:\n", fNonLinearityFunction);
2244 for(Int_t i=0; i<6; i++) printf("param[%d]=%f\n",i, fNonLinearityParams[i]);
2246 printf("Position Recalculation option %d, Particle Type %d, fW0 %2.2f, Recalibrate Data %d \n",fPosAlgo,fParticleType,fW0, fRecalibration);
2248 printf("Matching criteria: ");
2251 printf("sqrt(dEta^2+dPhi^2)<%4.3f\n",fCutR);
2253 else if(fCutEtaPhiSeparate)
2255 printf("dEta<%4.3f, dPhi<%4.3f\n",fCutEta,fCutPhi);
2260 printf("please specify your cut criteria\n");
2261 printf("To cut on sqrt(dEta^2+dPhi^2), use: SwitchOnCutEtaPhiSum()\n");
2262 printf("To cut on dEta and dPhi separately, use: SwitchOnCutEtaPhiSeparate()\n");
2265 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);
2266 printf("Cluster selection window: dR < %2.0f\n",fClusterWindow);
2268 printf("Track cuts: \n");
2269 printf("Minimum track pT: %1.2f\n",fCutMinTrackPt);
2270 printf("AOD track selection mask: %d\n",fAODFilterMask);
2271 printf("TPCRefit = %d, ITSRefit = %d\n",fCutRequireTPCRefit,fCutRequireITSRefit);
2272 printf("AcceptKinks = %d\n",fCutAcceptKinkDaughters);
2273 printf("MinNCulsterTPC = %d, MinNClusterITS = %d\n",fCutMinNClusterTPC,fCutMinNClusterITS);
2274 printf("MaxChi2TPC = %2.2f, MaxChi2ITS = %2.2f\n",fCutMaxChi2PerClusterTPC,fCutMaxChi2PerClusterITS);
2275 printf("DCSToVertex2D = %d, MaxDCAToVertexXY = %2.2f, MaxDCAToVertexZ = %2.2f\n",fCutDCAToVertex2D,fCutMaxDCAToVertexXY,fCutMaxDCAToVertexZ);
2278 //_____________________________________________________________________
2279 void AliEMCALRecoUtils::SetRunDependentCorrections(Int_t runnumber)
2281 //Get EMCAL time dependent corrections from file and put them in the recalibration histograms
2282 //Do it only once and only if it is requested
2284 if(!fUseRunCorrectionFactors) return;
2285 if(fRunCorrectionFactorsSet) return;
2287 AliInfo(Form("AliEMCALRecoUtils::GetRunDependentCorrections() - Get Correction Factors for Run number %d\n",runnumber));
2289 AliEMCALCalibTimeDepCorrection *corr = new AliEMCALCalibTimeDepCorrection();
2290 corr->ReadRootInfo(Form("CorrectionFiles/Run%d_Correction.root",runnumber));
2292 SwitchOnRecalibration();
2293 for(Int_t ism = 0; ism < 4; ism++){
2294 for(Int_t icol = 0; icol < 48; icol++){
2295 for(Int_t irow = 0; irow < 24; irow++){
2296 Float_t orgRecalFactor = GetEMCALChannelRecalibrationFactors(ism)->GetBinContent(icol,irow);
2297 Float_t newRecalFactor = orgRecalFactor*corr->GetCorrection(ism, icol,irow,0);
2298 GetEMCALChannelRecalibrationFactors(ism)->SetBinContent(icol,irow,newRecalFactor);
2299 //printf("ism %d, icol %d, irow %d, corrections : org %f, time dep %f, final %f (org*time %f)\n",ism, icol, irow,
2300 // orgRecalFactor, corr->GetCorrection(ism, icol,irow,0),
2301 // (GetEMCALChannelRecalibrationFactors(ism))->GetBinContent(icol,irow),newRecalFactor);
2305 fRunCorrectionFactorsSet = kTRUE;