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 "AliEMCALPIDUtils.h"
62 ClassImp(AliEMCALRecoUtils)
64 //_____________________________________
65 AliEMCALRecoUtils::AliEMCALRecoUtils():
66 fParticleType(0), fPosAlgo(0), fW0(0),
67 fNonLinearityFunction(0), fNonLinearThreshold(0),
68 fSmearClusterEnergy(kFALSE), fRandom(),
69 fCellsRecalibrated(kFALSE), fRecalibration(kFALSE), fEMCALRecalibrationFactors(),
70 fTimeRecalibration(kFALSE), fEMCALTimeRecalibrationFactors(), fUseRunCorrectionFactors(kFALSE),
71 fRemoveBadChannels(kFALSE), fRecalDistToBadChannels(kFALSE), fEMCALBadChannelMap(),
72 fNCellsFromEMCALBorder(0), fNoEMCALBorderAtEta0(kTRUE),
73 fRejectExoticCluster(kFALSE), fRejectExoticCells(kFALSE),
74 fExoticCellFraction(0), fExoticCellDiffTime(0), fExoticCellMinAmplitude(0),
75 fPIDUtils(), fAODFilterMask(0),
76 fMatchedTrackIndex(0x0), fMatchedClusterIndex(0x0),
77 fResidualEta(0x0), fResidualPhi(0x0), fCutEtaPhiSum(kFALSE), fCutEtaPhiSeparate(kFALSE),
78 fCutR(0), fCutEta(0), fCutPhi(0),
79 fClusterWindow(0), fMass(0),
80 fStepSurface(0), fStepCluster(0),
81 fTrackCutsType(0), fCutMinTrackPt(0), fCutMinNClusterTPC(0),
82 fCutMinNClusterITS(0), fCutMaxChi2PerClusterTPC(0), fCutMaxChi2PerClusterITS(0),
83 fCutRequireTPCRefit(kFALSE), fCutRequireITSRefit(kFALSE), fCutAcceptKinkDaughters(kFALSE),
84 fCutMaxDCAToVertexXY(0), fCutMaxDCAToVertexZ(0), fCutDCAToVertex2D(kFALSE)
88 // Initialize all constant values which have to be used
89 // during Reco algorithm execution
96 fMatchedTrackIndex = new TArrayI();
97 fMatchedClusterIndex = new TArrayI();
98 fResidualPhi = new TArrayF();
99 fResidualEta = new TArrayF();
100 fPIDUtils = new AliEMCALPIDUtils();
105 //______________________________________________________________________
106 AliEMCALRecoUtils::AliEMCALRecoUtils(const AliEMCALRecoUtils & reco)
108 fParticleType(reco.fParticleType), fPosAlgo(reco.fPosAlgo), fW0(reco.fW0),
109 fNonLinearityFunction(reco.fNonLinearityFunction), fNonLinearThreshold(reco.fNonLinearThreshold),
110 fSmearClusterEnergy(reco.fSmearClusterEnergy), fRandom(),
111 fCellsRecalibrated(reco.fCellsRecalibrated),
112 fRecalibration(reco.fRecalibration), fEMCALRecalibrationFactors(reco.fEMCALRecalibrationFactors),
113 fTimeRecalibration(reco.fTimeRecalibration), fEMCALTimeRecalibrationFactors(reco.fEMCALTimeRecalibrationFactors),
114 fUseRunCorrectionFactors(reco.fUseRunCorrectionFactors),
115 fRemoveBadChannels(reco.fRemoveBadChannels), fRecalDistToBadChannels(reco.fRecalDistToBadChannels),
116 fEMCALBadChannelMap(reco.fEMCALBadChannelMap),
117 fNCellsFromEMCALBorder(reco.fNCellsFromEMCALBorder), fNoEMCALBorderAtEta0(reco.fNoEMCALBorderAtEta0),
118 fRejectExoticCluster(reco.fRejectExoticCluster), fRejectExoticCells(reco.fRejectExoticCells),
119 fExoticCellFraction(reco.fExoticCellFraction), fExoticCellDiffTime(reco.fExoticCellDiffTime),
120 fExoticCellMinAmplitude(reco.fExoticCellMinAmplitude),
121 fPIDUtils(reco.fPIDUtils), fAODFilterMask(reco.fAODFilterMask),
122 fMatchedTrackIndex( reco.fMatchedTrackIndex? new TArrayI(*reco.fMatchedTrackIndex):0x0),
123 fMatchedClusterIndex(reco.fMatchedClusterIndex?new TArrayI(*reco.fMatchedClusterIndex):0x0),
124 fResidualEta( reco.fResidualEta? new TArrayF(*reco.fResidualEta):0x0),
125 fResidualPhi( reco.fResidualPhi? new TArrayF(*reco.fResidualPhi):0x0),
126 fCutEtaPhiSum(reco.fCutEtaPhiSum), fCutEtaPhiSeparate(reco.fCutEtaPhiSeparate),
127 fCutR(reco.fCutR), fCutEta(reco.fCutEta), fCutPhi(reco.fCutPhi),
128 fClusterWindow(reco.fClusterWindow),
129 fMass(reco.fMass), fStepSurface(reco.fStepSurface), fStepCluster(reco.fStepCluster),
130 fTrackCutsType(reco.fTrackCutsType), fCutMinTrackPt(reco.fCutMinTrackPt),
131 fCutMinNClusterTPC(reco.fCutMinNClusterTPC), fCutMinNClusterITS(reco.fCutMinNClusterITS),
132 fCutMaxChi2PerClusterTPC(reco.fCutMaxChi2PerClusterTPC), fCutMaxChi2PerClusterITS(reco.fCutMaxChi2PerClusterITS),
133 fCutRequireTPCRefit(reco.fCutRequireTPCRefit), fCutRequireITSRefit(reco.fCutRequireITSRefit),
134 fCutAcceptKinkDaughters(reco.fCutAcceptKinkDaughters), fCutMaxDCAToVertexXY(reco.fCutMaxDCAToVertexXY),
135 fCutMaxDCAToVertexZ(reco.fCutMaxDCAToVertexZ), fCutDCAToVertex2D(reco.fCutDCAToVertex2D)
139 for(Int_t i = 0; i < 15 ; i++) { fMisalRotShift[i] = reco.fMisalRotShift[i] ;
140 fMisalTransShift[i] = reco.fMisalTransShift[i] ; }
141 for(Int_t i = 0; i < 7 ; i++) { fNonLinearityParams[i] = reco.fNonLinearityParams[i] ; }
142 for(Int_t i = 0; i < 3 ; i++) { fSmearClusterParam[i] = reco.fSmearClusterParam[i] ; }
147 //______________________________________________________________________
148 AliEMCALRecoUtils & AliEMCALRecoUtils::operator = (const AliEMCALRecoUtils & reco)
150 //Assignment operator
152 if(this == &reco)return *this;
153 ((TNamed *)this)->operator=(reco);
155 for(Int_t i = 0; i < 15 ; i++) { fMisalTransShift[i] = reco.fMisalTransShift[i] ;
156 fMisalRotShift[i] = reco.fMisalRotShift[i] ; }
157 for(Int_t i = 0; i < 7 ; i++) { fNonLinearityParams[i] = reco.fNonLinearityParams[i] ; }
158 for(Int_t i = 0; i < 3 ; i++) { fSmearClusterParam[i] = reco.fSmearClusterParam[i] ; }
160 fParticleType = reco.fParticleType;
161 fPosAlgo = reco.fPosAlgo;
164 fNonLinearityFunction = reco.fNonLinearityFunction;
165 fNonLinearThreshold = reco.fNonLinearThreshold;
166 fSmearClusterEnergy = reco.fSmearClusterEnergy;
168 fCellsRecalibrated = reco.fCellsRecalibrated;
169 fRecalibration = reco.fRecalibration;
170 fEMCALRecalibrationFactors = reco.fEMCALRecalibrationFactors;
172 fTimeRecalibration = reco.fTimeRecalibration;
173 fEMCALTimeRecalibrationFactors = reco.fEMCALTimeRecalibrationFactors;
175 fUseRunCorrectionFactors = reco.fUseRunCorrectionFactors;
177 fRemoveBadChannels = reco.fRemoveBadChannels;
178 fRecalDistToBadChannels = reco.fRecalDistToBadChannels;
179 fEMCALBadChannelMap = reco.fEMCALBadChannelMap;
181 fNCellsFromEMCALBorder = reco.fNCellsFromEMCALBorder;
182 fNoEMCALBorderAtEta0 = reco.fNoEMCALBorderAtEta0;
184 fRejectExoticCluster = reco.fRejectExoticCluster;
185 fRejectExoticCells = reco.fRejectExoticCells;
186 fExoticCellFraction = reco.fExoticCellFraction;
187 fExoticCellDiffTime = reco.fExoticCellDiffTime;
188 fExoticCellMinAmplitude = reco.fExoticCellMinAmplitude;
190 fPIDUtils = reco.fPIDUtils;
192 fAODFilterMask = reco.fAODFilterMask;
194 fCutEtaPhiSum = reco.fCutEtaPhiSum;
195 fCutEtaPhiSeparate = reco.fCutEtaPhiSeparate;
197 fCutEta = reco.fCutEta;
198 fCutPhi = reco.fCutPhi;
199 fClusterWindow = reco.fClusterWindow;
201 fStepSurface = reco.fStepSurface;
202 fStepCluster = reco.fStepCluster;
204 fTrackCutsType = reco.fTrackCutsType;
205 fCutMinTrackPt = reco.fCutMinTrackPt;
206 fCutMinNClusterTPC = reco.fCutMinNClusterTPC;
207 fCutMinNClusterITS = reco.fCutMinNClusterITS;
208 fCutMaxChi2PerClusterTPC = reco.fCutMaxChi2PerClusterTPC;
209 fCutMaxChi2PerClusterITS = reco.fCutMaxChi2PerClusterITS;
210 fCutRequireTPCRefit = reco.fCutRequireTPCRefit;
211 fCutRequireITSRefit = reco.fCutRequireITSRefit;
212 fCutAcceptKinkDaughters = reco.fCutAcceptKinkDaughters;
213 fCutMaxDCAToVertexXY = reco.fCutMaxDCAToVertexXY;
214 fCutMaxDCAToVertexZ = reco.fCutMaxDCAToVertexZ;
215 fCutDCAToVertex2D = reco.fCutDCAToVertex2D;
217 if(reco.fResidualEta)
219 // assign or copy construct
222 *fResidualEta = *reco.fResidualEta;
226 fResidualEta = new TArrayF(*reco.fResidualEta);
231 if(fResidualEta)delete fResidualEta;
235 if(reco.fResidualPhi)
237 // assign or copy construct
240 *fResidualPhi = *reco.fResidualPhi;
244 fResidualPhi = new TArrayF(*reco.fResidualPhi);
249 if(fResidualPhi)delete fResidualPhi;
253 if(reco.fMatchedTrackIndex)
255 // assign or copy construct
256 if(fMatchedTrackIndex)
258 *fMatchedTrackIndex = *reco.fMatchedTrackIndex;
262 fMatchedTrackIndex = new TArrayI(*reco.fMatchedTrackIndex);
267 if(fMatchedTrackIndex)delete fMatchedTrackIndex;
268 fMatchedTrackIndex = 0;
271 if(reco.fMatchedClusterIndex)
273 // assign or copy construct
274 if(fMatchedClusterIndex)
276 *fMatchedClusterIndex = *reco.fMatchedClusterIndex;
280 fMatchedClusterIndex = new TArrayI(*reco.fMatchedClusterIndex);
285 if(fMatchedClusterIndex)delete fMatchedClusterIndex;
286 fMatchedClusterIndex = 0;
293 //_____________________________________
294 AliEMCALRecoUtils::~AliEMCALRecoUtils()
298 if(fEMCALRecalibrationFactors)
300 fEMCALRecalibrationFactors->Clear();
301 delete fEMCALRecalibrationFactors;
304 if(fEMCALTimeRecalibrationFactors)
306 fEMCALTimeRecalibrationFactors->Clear();
307 delete fEMCALTimeRecalibrationFactors;
310 if(fEMCALBadChannelMap)
312 fEMCALBadChannelMap->Clear();
313 delete fEMCALBadChannelMap;
316 delete fMatchedTrackIndex ;
317 delete fMatchedClusterIndex ;
318 delete fResidualEta ;
319 delete fResidualPhi ;
325 //_______________________________________________________________________________
326 Bool_t AliEMCALRecoUtils::AcceptCalibrateCell(const Int_t absID, const Int_t bc,
327 Float_t & amp, Double_t & time,
328 AliVCaloCells* cells)
330 // Reject cell if criteria not passed and calibrate it
332 AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
334 if(absID < 0 || absID >= 24*48*geom->GetNumberOfSuperModules()) return kFALSE;
336 Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1;
338 if(!geom->GetCellIndex(absID,imod,iTower,iIphi,iIeta))
340 // cell absID does not exist
345 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);
347 // Do not include bad channels found in analysis,
348 if( IsBadChannelsRemovalSwitchedOn() && GetEMCALChannelStatus(imod, ieta, iphi))
354 amp = cells->GetCellAmplitude(absID);
355 if(!fCellsRecalibrated && IsRecalibrationOn())
356 amp *= GetEMCALChannelRecalibrationFactor(imod,ieta,iphi);
360 time = cells->GetCellTime(absID);
362 RecalibrateCellTime(absID,bc,time);
367 //_____________________________________________________________________________
368 Bool_t AliEMCALRecoUtils::CheckCellFiducialRegion(const AliEMCALGeometry* geom,
369 const AliVCluster* cluster,
370 AliVCaloCells* cells)
372 // Given the list of AbsId of the cluster, get the maximum cell and
373 // check if there are fNCellsFromBorder from the calorimeter border
377 AliInfo("Cluster pointer null!");
381 //If the distance to the border is 0 or negative just exit accept all clusters
382 if(cells->GetType()==AliVCaloCells::kEMCALCell && fNCellsFromEMCALBorder <= 0 ) return kTRUE;
384 Int_t absIdMax = -1, iSM =-1, ieta = -1, iphi = -1;
385 Bool_t shared = kFALSE;
386 GetMaxEnergyCell(geom, cells, cluster, absIdMax, iSM, ieta, iphi, shared);
388 AliDebug(2,Form("Cluster Max AbsId %d, Cell Energy %2.2f, Cluster Energy %2.2f, Ncells from border %d, EMCAL eta=0 %d\n",
389 absIdMax, cells->GetCellAmplitude(absIdMax), cluster->E(), fNCellsFromEMCALBorder, fNoEMCALBorderAtEta0));
391 if(absIdMax==-1) return kFALSE;
393 //Check if the cell is close to the borders:
394 Bool_t okrow = kFALSE;
395 Bool_t okcol = kFALSE;
397 if(iSM < 0 || iphi < 0 || ieta < 0 )
399 AliFatal(Form("Negative value for super module: %d, or cell ieta: %d, or cell iphi: %d, check EMCAL geometry name\n",
406 if(iphi >= fNCellsFromEMCALBorder && iphi < 24-fNCellsFromEMCALBorder) okrow =kTRUE;
408 else if (iSM >=10 && ( ( geom->GetEMCGeometry()->GetGeoName()).Contains("12SMV1")))
410 if(iphi >= fNCellsFromEMCALBorder && iphi < 8-fNCellsFromEMCALBorder) okrow =kTRUE; //1/3 sm case
414 if(iphi >= fNCellsFromEMCALBorder && iphi < 12-fNCellsFromEMCALBorder) okrow =kTRUE; // half SM case
418 if(!fNoEMCALBorderAtEta0)
420 if(ieta > fNCellsFromEMCALBorder && ieta < 48-fNCellsFromEMCALBorder) okcol =kTRUE;
426 if(ieta >= fNCellsFromEMCALBorder) okcol = kTRUE;
430 if(ieta < 48-fNCellsFromEMCALBorder) okcol = kTRUE;
434 AliDebug(2,Form("EMCAL Cluster in %d cells fiducial volume: ieta %d, iphi %d, SM %d: column? %d, row? %d\nq",
435 fNCellsFromEMCALBorder, ieta, iphi, iSM, okcol, okrow));
439 //printf("Accept\n");
444 //printf("Reject\n");
445 AliDebug(2,Form("Reject cluster in border, max cell : ieta %d, iphi %d, SM %d\n",ieta, iphi, iSM));
452 //_______________________________________________________________________________
453 Bool_t AliEMCALRecoUtils::ClusterContainsBadChannel(const AliEMCALGeometry* geom,
454 const UShort_t* cellList,
457 // Check that in the cluster cells, there is no bad channel of those stored
458 // in fEMCALBadChannelMap or fPHOSBadChannelMap
460 if(!fRemoveBadChannels) return kFALSE;
461 if(!fEMCALBadChannelMap) return kFALSE;
466 for(Int_t iCell = 0; iCell<nCells; iCell++)
468 //Get the column and row
469 Int_t iTower = -1, iIphi = -1, iIeta = -1;
470 geom->GetCellIndex(cellList[iCell],imod,iTower,iIphi,iIeta);
471 if(fEMCALBadChannelMap->GetEntries() <= imod) continue;
472 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,irow,icol);
473 if(GetEMCALChannelStatus(imod, icol, irow))
475 AliDebug(2,Form("Cluster with bad channel: SM %d, col %d, row %d\n",imod, icol, irow));
479 }// cell cluster loop
485 //___________________________________________________________________________
486 Float_t AliEMCALRecoUtils::GetECross(const Int_t absID, const Double_t tcell,
487 AliVCaloCells* cells, const Int_t bc)
489 //Calculate the energy in the cross around the energy given cell
491 AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
493 Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1;
494 geom->GetCellIndex(absID,imod,iTower,iIphi,iIeta);
495 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);
497 //Get close cells index, energy and time, not in corners
502 if( iphi < AliEMCALGeoParams::fgkEMCALRows-1) absID1 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi+1, ieta);
503 if( iphi > 0 ) absID2 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi-1, ieta);
505 // 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+1, 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-1, iphi, AliEMCALGeoParams::fgkEMCALCols-1);
522 if( ieta < AliEMCALGeoParams::fgkEMCALCols-1 )
523 absID3 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta+1);
525 absID4 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta-1);
528 //printf("IMOD %d, AbsId %d, a %d, b %d, c %d e %d \n",imod,absID,absID1,absID2,absID3,absID4);
530 Float_t ecell1 = 0, ecell2 = 0, ecell3 = 0, ecell4 = 0;
531 Double_t tcell1 = 0, tcell2 = 0, tcell3 = 0, tcell4 = 0;
532 Bool_t accept1 = 0, accept2 = 0, accept3 = 0, accept4 = 0;
534 accept1 = AcceptCalibrateCell(absID1,bc, ecell1,tcell1,cells);
535 accept2 = AcceptCalibrateCell(absID2,bc, ecell2,tcell2,cells);
536 accept3 = AcceptCalibrateCell(absID3,bc, ecell3,tcell3,cells);
537 accept4 = AcceptCalibrateCell(absID4,bc, ecell4,tcell4,cells);
540 printf("Cell absID %d \n",absID);
541 printf("\t accept1 %d, accept2 %d, accept3 %d, accept4 %d\n",
542 accept1,accept2,accept3,accept4);
543 printf("\t id %d: id1 %d, id2 %d, id3 %d, id4 %d\n",
544 absID,absID1,absID2,absID3,absID4);
545 printf("\t e %f: e1 %f, e2 %f, e3 %f, e4 %f\n",
546 ecell,ecell1,ecell2,ecell3,ecell4);
547 printf("\t t %f: t1 %f, t2 %f, t3 %f, t4 %f;\n dt1 %f, dt2 %f, dt3 %f, dt4 %f\n",
548 tcell*1.e9,tcell1*1.e9,tcell2*1.e9,tcell3*1.e9,tcell4*1.e9,
549 TMath::Abs(tcell-tcell1)*1.e9, TMath::Abs(tcell-tcell2)*1.e9, TMath::Abs(tcell-tcell3)*1.e9, TMath::Abs(tcell-tcell4)*1.e9);
552 if(TMath::Abs(tcell-tcell1)*1.e9 > fExoticCellDiffTime) ecell1 = 0 ;
553 if(TMath::Abs(tcell-tcell2)*1.e9 > fExoticCellDiffTime) ecell2 = 0 ;
554 if(TMath::Abs(tcell-tcell3)*1.e9 > fExoticCellDiffTime) ecell3 = 0 ;
555 if(TMath::Abs(tcell-tcell4)*1.e9 > fExoticCellDiffTime) ecell4 = 0 ;
557 return ecell1+ecell2+ecell3+ecell4;
561 //_____________________________________________________________________________________________
562 Bool_t AliEMCALRecoUtils::IsExoticCell(const Int_t absID, AliVCaloCells* cells, const Int_t bc)
564 // Look to cell neighbourhood and reject if it seems exotic
565 // Do before recalibrating the cells
567 if(!fRejectExoticCells) return kFALSE;
571 Bool_t accept = AcceptCalibrateCell(absID, bc, ecell ,tcell ,cells);
573 if(!accept) return kTRUE; // reject this cell
575 if(ecell < fExoticCellMinAmplitude) return kFALSE; // do not reject low energy cells
577 Float_t eCross = GetECross(absID,tcell,cells,bc);
579 if(1-eCross/ecell > fExoticCellFraction)
581 AliDebug(2,Form("AliEMCALRecoUtils::IsExoticCell() - EXOTIC CELL id %d, eCell %f, eCross %f, 1-eCross/eCell %f\n",
582 absID,ecell,eCross,1-eCross/ecell));
589 //___________________________________________________________________
590 Bool_t AliEMCALRecoUtils::IsExoticCluster(const AliVCluster *cluster,
591 AliVCaloCells *cells,
594 // Check if the cluster highest energy tower is exotic
598 AliInfo("Cluster pointer null!");
602 if(!fRejectExoticCluster) return kFALSE;
604 // Get highest energy tower
605 AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
606 Int_t iSupMod = -1, absId = -1, ieta = -1, iphi = -1;
607 Bool_t shared = kFALSE;
608 GetMaxEnergyCell(geom, cells, cluster, absId, iSupMod, ieta, iphi, shared);
610 return IsExoticCell(absId,cells,bc);
614 //_______________________________________________________________________
615 Float_t AliEMCALRecoUtils::SmearClusterEnergy(const AliVCluster* cluster)
617 //In case of MC analysis, smear energy to match resolution/calibration in real data
621 AliInfo("Cluster pointer null!");
625 Float_t energy = cluster->E() ;
626 Float_t rdmEnergy = energy ;
627 if(fSmearClusterEnergy)
629 rdmEnergy = fRandom.Gaus(energy,fSmearClusterParam[0] * TMath::Sqrt(energy) +
630 fSmearClusterParam[1] * energy +
631 fSmearClusterParam[2] );
632 AliDebug(2, Form("Energy: original %f, smeared %f\n", energy, rdmEnergy));
638 //____________________________________________________________________________
639 Float_t AliEMCALRecoUtils::CorrectClusterEnergyLinearity(AliVCluster* cluster)
641 // Correct cluster energy from non linearity functions
645 AliInfo("Cluster pointer null!");
649 Float_t energy = cluster->E();
651 switch (fNonLinearityFunction)
656 //Non-Linearity correction (from MC with function ([0]*exp(-[1]/E))+(([2]/([3]*2.*TMath::Pi())*exp(-(E-[4])^2/(2.*[3]^2)))))
657 //fNonLinearityParams[0] = 1.014;
658 //fNonLinearityParams[1] =-0.03329;
659 //fNonLinearityParams[2] =-0.3853;
660 //fNonLinearityParams[3] = 0.5423;
661 //fNonLinearityParams[4] =-0.4335;
662 energy *= (fNonLinearityParams[0]*exp(-fNonLinearityParams[1]/energy))+
663 ((fNonLinearityParams[2]/(fNonLinearityParams[3]*2.*TMath::Pi())*
664 exp(-(energy-fNonLinearityParams[4])*(energy-fNonLinearityParams[4])/(2.*fNonLinearityParams[3]*fNonLinearityParams[3]))));
670 //Non-Linearity correction (from MC with function [0]/((x+[1])^[2]))+1;
671 //fNonLinearityParams[0] = 3.11111e-02;
672 //fNonLinearityParams[1] =-5.71666e-02;
673 //fNonLinearityParams[2] = 5.67995e-01;
675 energy *= fNonLinearityParams[0]/TMath::Power(energy+fNonLinearityParams[1],fNonLinearityParams[2])+1;
681 //Non-Linearity correction (from Olga Data with function p0+p1*exp(-p2*E))
682 //fNonLinearityParams[0] = 1.04;
683 //fNonLinearityParams[1] = -0.1445;
684 //fNonLinearityParams[2] = 1.046;
685 energy /= (fNonLinearityParams[0]+fNonLinearityParams[1]*exp(-fNonLinearityParams[2]*energy)); //Olga function
689 case kPi0GammaConversion:
691 //Non-Linearity correction (Nicolas from Dimitri Data with function C*[1-a*exp(-b*E)])
692 //fNonLinearityParams[0] = 0.139393/0.1349766;
693 //fNonLinearityParams[1] = 0.0566186;
694 //fNonLinearityParams[2] = 0.982133;
695 energy /= fNonLinearityParams[0]*(1-fNonLinearityParams[1]*exp(-fNonLinearityParams[2]*energy));
702 //From beam test, Alexei's results, for different ZS thresholds
703 // th=30 MeV; th = 45 MeV; th = 75 MeV
704 //fNonLinearityParams[0] = 1.007; 1.003; 1.002
705 //fNonLinearityParams[1] = 0.894; 0.719; 0.797
706 //fNonLinearityParams[2] = 0.246; 0.334; 0.358
707 //Rescale the param[0] with 1.03
708 energy /= fNonLinearityParams[0]/(1+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]));
713 case kBeamTestCorrected:
715 //From beam test, corrected for material between beam and EMCAL
716 //fNonLinearityParams[0] = 0.99078
717 //fNonLinearityParams[1] = 0.161499;
718 //fNonLinearityParams[2] = 0.655166;
719 //fNonLinearityParams[3] = 0.134101;
720 //fNonLinearityParams[4] = 163.282;
721 //fNonLinearityParams[5] = 23.6904;
722 //fNonLinearityParams[6] = 0.978;
723 energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
729 AliDebug(2,"No correction on the energy\n");
737 //__________________________________________________
738 void AliEMCALRecoUtils::InitNonLinearityParam()
740 //Initialising Non Linearity Parameters
742 if(fNonLinearityFunction == kPi0MC)
744 fNonLinearityParams[0] = 1.014;
745 fNonLinearityParams[1] = -0.03329;
746 fNonLinearityParams[2] = -0.3853;
747 fNonLinearityParams[3] = 0.5423;
748 fNonLinearityParams[4] = -0.4335;
751 if(fNonLinearityFunction == kPi0GammaGamma)
753 fNonLinearityParams[0] = 1.04;
754 fNonLinearityParams[1] = -0.1445;
755 fNonLinearityParams[2] = 1.046;
758 if(fNonLinearityFunction == kPi0GammaConversion)
760 fNonLinearityParams[0] = 0.139393;
761 fNonLinearityParams[1] = 0.0566186;
762 fNonLinearityParams[2] = 0.982133;
765 if(fNonLinearityFunction == kBeamTest)
767 if(fNonLinearThreshold == 30)
769 fNonLinearityParams[0] = 1.007;
770 fNonLinearityParams[1] = 0.894;
771 fNonLinearityParams[2] = 0.246;
773 if(fNonLinearThreshold == 45)
775 fNonLinearityParams[0] = 1.003;
776 fNonLinearityParams[1] = 0.719;
777 fNonLinearityParams[2] = 0.334;
779 if(fNonLinearThreshold == 75)
781 fNonLinearityParams[0] = 1.002;
782 fNonLinearityParams[1] = 0.797;
783 fNonLinearityParams[2] = 0.358;
787 if(fNonLinearityFunction == kBeamTestCorrected)
789 fNonLinearityParams[0] = 0.99078;
790 fNonLinearityParams[1] = 0.161499;
791 fNonLinearityParams[2] = 0.655166;
792 fNonLinearityParams[3] = 0.134101;
793 fNonLinearityParams[4] = 163.282;
794 fNonLinearityParams[5] = 23.6904;
795 fNonLinearityParams[6] = 0.978;
799 //_________________________________________________________
800 Float_t AliEMCALRecoUtils::GetDepth(const Float_t energy,
801 const Int_t iParticle,
802 const Int_t iSM) const
804 //Calculate shower depth for a given cluster energy and particle type
814 depth = x0 * (TMath::Log(energy*1000/ ecr) + 0.5); //Multiply energy by 1000 to transform to MeV
818 depth = x0 * (TMath::Log(energy*1000/ ecr) - 0.5); //Multiply energy by 1000 to transform to MeV
826 gGeoManager->cd("ALIC_1/XEN1_1");
827 TGeoNode *geoXEn1 = gGeoManager->GetCurrentNode();
828 TGeoNodeMatrix *geoSM = dynamic_cast<TGeoNodeMatrix *>(geoXEn1->GetDaughter(iSM));
831 TGeoVolume *geoSMVol = geoSM->GetVolume();
832 TGeoShape *geoSMShape = geoSMVol->GetShape();
833 TGeoBBox *geoBox = dynamic_cast<TGeoBBox *>(geoSMShape);
834 if(geoBox) depth = 0.5 * geoBox->GetDX()*2 ;
835 else AliFatal("Null GEANT box");
837 else AliFatal("NULL GEANT node matrix");
841 depth = x0 * (TMath::Log(energy*1000 / ecr) - 0.5); //Multiply energy by 1000 to transform to MeV
847 depth = x0 * (TMath::Log(energy*1000 / ecr) + 0.5); //Multiply energy by 1000 to transform to MeV
853 //____________________________________________________________________
854 void AliEMCALRecoUtils::GetMaxEnergyCell(const AliEMCALGeometry *geom,
855 AliVCaloCells* cells,
856 const AliVCluster* clu,
863 //For a given CaloCluster gets the absId of the cell
864 //with maximum energy deposit.
867 Double_t eCell = -1.;
868 Float_t fraction = 1.;
869 Float_t recalFactor = 1.;
870 Int_t cellAbsId = -1 ;
879 AliInfo("Cluster pointer null!");
880 absId=-1; iSupMod0=-1, ieta = -1; iphi = -1; shared = -1;
884 for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++)
886 cellAbsId = clu->GetCellAbsId(iDig);
887 fraction = clu->GetCellAmplitudeFraction(iDig);
888 //printf("a Cell %d, id, %d, amp %f, fraction %f\n",iDig,cellAbsId,cells->GetCellAmplitude(cellAbsId),fraction);
889 if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
890 geom->GetCellIndex(cellAbsId,iSupMod,iTower,iIphi,iIeta);
891 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);
896 else if(iSupMod0!=iSupMod)
899 //printf("AliEMCALRecoUtils::GetMaxEnergyCell() - SHARED CLUSTER\n");
901 if(!fCellsRecalibrated && IsRecalibrationOn())
903 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
905 eCell = cells->GetCellAmplitude(cellAbsId)*fraction*recalFactor;
906 //printf("b Cell %d, id, %d, amp %f, fraction %f\n",iDig,cellAbsId,eCell,fraction);
911 //printf("\t new max: cell %d, e %f, ecell %f\n",maxId, eMax,eCell);
915 //Get from the absid the supermodule, tower and eta/phi numbers
916 geom->GetCellIndex(absId,iSupMod,iTower,iIphi,iIeta);
917 //Gives SuperModule and Tower numbers
918 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,
919 iIphi, iIeta,iphi,ieta);
920 //printf("Max id %d, iSM %d, col %d, row %d\n",absId,iSupMod,ieta,iphi);
921 //printf("Max end---\n");
924 //______________________________________
925 void AliEMCALRecoUtils::InitParameters()
927 // Initialize data members with default values
929 fParticleType = kPhoton;
930 fPosAlgo = kUnchanged;
933 fNonLinearityFunction = kNoCorrection;
934 fNonLinearThreshold = 30;
936 fExoticCellFraction = 0.97;
937 fExoticCellDiffTime = 1e6;
938 fExoticCellMinAmplitude = 0.5;
942 fCutEtaPhiSum = kTRUE;
943 fCutEtaPhiSeparate = kFALSE;
949 fClusterWindow = 100;
954 fTrackCutsType = kLooseCut;
957 fCutMinNClusterTPC = -1;
958 fCutMinNClusterITS = -1;
960 fCutMaxChi2PerClusterTPC = 1e10;
961 fCutMaxChi2PerClusterITS = 1e10;
963 fCutRequireTPCRefit = kFALSE;
964 fCutRequireITSRefit = kFALSE;
965 fCutAcceptKinkDaughters = kFALSE;
967 fCutMaxDCAToVertexXY = 1e10;
968 fCutMaxDCAToVertexZ = 1e10;
969 fCutDCAToVertex2D = kFALSE;
972 //Misalignment matrices
973 for(Int_t i = 0; i < 15 ; i++)
975 fMisalTransShift[i] = 0.;
976 fMisalRotShift[i] = 0.;
980 for(Int_t i = 0; i < 7 ; i++) fNonLinearityParams[i] = 0.;
982 //For kBeamTestCorrected case, but default is no correction
983 fNonLinearityParams[0] = 0.99078;
984 fNonLinearityParams[1] = 0.161499;
985 fNonLinearityParams[2] = 0.655166;
986 fNonLinearityParams[3] = 0.134101;
987 fNonLinearityParams[4] = 163.282;
988 fNonLinearityParams[5] = 23.6904;
989 fNonLinearityParams[6] = 0.978;
991 //For kPi0GammaGamma case
992 //fNonLinearityParams[0] = 0.1457/0.1349766/1.038;
993 //fNonLinearityParams[1] = -0.02024/0.1349766/1.038;
994 //fNonLinearityParams[2] = 1.046;
996 //Cluster energy smearing
997 fSmearClusterEnergy = kFALSE;
998 fSmearClusterParam[0] = 0.07; // * sqrt E term
999 fSmearClusterParam[1] = 0.00; // * E term
1000 fSmearClusterParam[2] = 0.00; // constant
1003 //_____________________________________________________
1004 void AliEMCALRecoUtils::InitEMCALRecalibrationFactors()
1006 //Init EMCAL recalibration factors
1007 AliDebug(2,"AliCalorimeterUtils::InitEMCALRecalibrationFactors()");
1008 //In order to avoid rewriting the same histograms
1009 Bool_t oldStatus = TH1::AddDirectoryStatus();
1010 TH1::AddDirectory(kFALSE);
1012 fEMCALRecalibrationFactors = new TObjArray(12);
1013 for (int i = 0; i < 12; i++)
1014 fEMCALRecalibrationFactors->Add(new TH2F(Form("EMCALRecalFactors_SM%d",i),
1015 Form("EMCALRecalFactors_SM%d",i), 48, 0, 48, 24, 0, 24));
1016 //Init the histograms with 1
1017 for (Int_t sm = 0; sm < 12; sm++)
1019 for (Int_t i = 0; i < 48; i++)
1021 for (Int_t j = 0; j < 24; j++)
1023 SetEMCALChannelRecalibrationFactor(sm,i,j,1.);
1028 fEMCALRecalibrationFactors->SetOwner(kTRUE);
1029 fEMCALRecalibrationFactors->Compress();
1031 //In order to avoid rewriting the same histograms
1032 TH1::AddDirectory(oldStatus);
1035 //_________________________________________________________
1036 void AliEMCALRecoUtils::InitEMCALTimeRecalibrationFactors()
1038 //Init EMCAL recalibration factors
1039 AliDebug(2,"AliCalorimeterUtils::InitEMCALRecalibrationFactors()");
1040 //In order to avoid rewriting the same histograms
1041 Bool_t oldStatus = TH1::AddDirectoryStatus();
1042 TH1::AddDirectory(kFALSE);
1044 fEMCALTimeRecalibrationFactors = new TObjArray(4);
1045 for (int i = 0; i < 4; i++)
1046 fEMCALTimeRecalibrationFactors->Add(new TH1F(Form("hAllTimeAvBC%d",i),
1047 Form("hAllTimeAvBC%d",i),
1048 48*24*12,0.,48*24*12) );
1049 //Init the histograms with 1
1050 for (Int_t bc = 0; bc < 4; bc++)
1052 for (Int_t i = 0; i < 48*24*12; i++)
1053 SetEMCALChannelTimeRecalibrationFactor(bc,i,0.);
1056 fEMCALTimeRecalibrationFactors->SetOwner(kTRUE);
1057 fEMCALTimeRecalibrationFactors->Compress();
1059 //In order to avoid rewriting the same histograms
1060 TH1::AddDirectory(oldStatus);
1063 //____________________________________________________
1064 void AliEMCALRecoUtils::InitEMCALBadChannelStatusMap()
1066 //Init EMCAL bad channels map
1067 AliDebug(2,"AliEMCALRecoUtils::InitEMCALBadChannelStatusMap()");
1068 //In order to avoid rewriting the same histograms
1069 Bool_t oldStatus = TH1::AddDirectoryStatus();
1070 TH1::AddDirectory(kFALSE);
1072 fEMCALBadChannelMap = new TObjArray(12);
1073 //TH2F * hTemp = new TH2I("EMCALBadChannelMap","EMCAL SuperModule bad channel map", 48, 0, 48, 24, 0, 24);
1074 for (int i = 0; i < 12; i++)
1076 fEMCALBadChannelMap->Add(new TH2I(Form("EMCALBadChannelMap_Mod%d",i),Form("EMCALBadChannelMap_Mod%d",i), 48, 0, 48, 24, 0, 24));
1079 fEMCALBadChannelMap->SetOwner(kTRUE);
1080 fEMCALBadChannelMap->Compress();
1082 //In order to avoid rewriting the same histograms
1083 TH1::AddDirectory(oldStatus);
1086 //____________________________________________________________________________
1087 void AliEMCALRecoUtils::RecalibrateClusterEnergy(const AliEMCALGeometry* geom,
1088 AliVCluster * cluster,
1089 AliVCaloCells * cells,
1092 // Recalibrate the cluster energy and Time, considering the recalibration map
1093 // and the energy of the cells and time that compose the cluster.
1094 // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1098 AliInfo("Cluster pointer null!");
1102 //Get the cluster number of cells and list of absId, check what kind of cluster do we have.
1103 UShort_t * index = cluster->GetCellsAbsId() ;
1104 Double_t * fraction = cluster->GetCellsAmplitudeFraction() ;
1105 Int_t ncells = cluster->GetNCells();
1107 //Initialize some used variables
1110 Int_t icol =-1, irow =-1, imod=1;
1111 Float_t factor = 1, frac = 0;
1112 Int_t absIdMax = -1;
1115 //Loop on the cells, get the cell amplitude and recalibration factor, multiply and and to the new energy
1116 for(Int_t icell = 0; icell < ncells; icell++)
1118 absId = index[icell];
1119 frac = fraction[icell];
1120 if(frac < 1e-5) frac = 1; //in case of EMCAL, this is set as 0 since unfolding is off
1122 if(!fCellsRecalibrated && IsRecalibrationOn())
1125 Int_t iTower = -1, iIphi = -1, iIeta = -1;
1126 geom->GetCellIndex(absId,imod,iTower,iIphi,iIeta);
1127 if(fEMCALRecalibrationFactors->GetEntries() <= imod) continue;
1128 geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,irow,icol);
1129 factor = GetEMCALChannelRecalibrationFactor(imod,icol,irow);
1131 AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - recalibrate cell: module %d, col %d, row %d, cell fraction %f,recalibration factor %f, cell energy %f\n",
1132 imod,icol,irow,frac,factor,cells->GetCellAmplitude(absId)));
1136 energy += cells->GetCellAmplitude(absId)*factor*frac;
1138 if(emax < cells->GetCellAmplitude(absId)*factor*frac)
1140 emax = cells->GetCellAmplitude(absId)*factor*frac;
1145 AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - Energy before %f, after %f \n",cluster->E(),energy));
1147 cluster->SetE(energy);
1149 // Recalculate time of cluster
1150 Double_t timeorg = cluster->GetTOF();
1151 if(!fCellsRecalibrated && IsTimeRecalibrationOn())
1153 Double_t time = timeorg;
1154 RecalibrateCellTime(absIdMax,bc,time);
1155 cluster->SetTOF(time);
1158 AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - Time before %f, after %f \n",timeorg,cluster->GetTOF()));
1162 //_____________________________________________________________
1163 void AliEMCALRecoUtils::RecalibrateCells(AliVCaloCells * cells,
1166 // Recalibrate the cells time and energy, considering the recalibration map and the energy
1167 // of the cells that compose the cluster.
1168 // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1170 if(!IsRecalibrationOn() && !IsTimeRecalibrationOn()) return;
1174 AliInfo("Cells pointer null!");
1179 Bool_t accept = kFALSE;
1182 Double_t ecellin = 0;
1183 Double_t tcellin = 0;
1184 Short_t mclabel = -1;
1187 Int_t nEMcell = cells->GetNumberOfCells() ;
1188 for (Int_t iCell = 0; iCell < nEMcell; iCell++)
1190 cells->GetCell( iCell, absId, ecellin, tcellin, mclabel, efrac );
1192 accept = AcceptCalibrateCell(absId, bc, ecell ,tcell ,cells);
1200 cells->SetCell(iCell,absId,ecell, tcell, mclabel, efrac);
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,
1507 Float_t & l0, Float_t & l1,
1508 Float_t & disp, Float_t & dEta, Float_t & dPhi,
1509 Float_t & sEta, Float_t & sPhi, Float_t & sEtaPhi)
1511 // Calculates new center of gravity in the local EMCAL-module coordinates
1512 // and tranfers into global ALICE coordinates
1513 // Calculates Dispersion and main axis
1517 AliInfo("Cluster pointer null!");
1521 Double_t eCell = 0.;
1522 Float_t fraction = 1.;
1523 Float_t recalFactor = 1.;
1531 Double_t etai = -1.;
1532 Double_t phii = -1.;
1537 Double_t etaMean = 0.;
1538 Double_t phiMean = 0.;
1541 for(Int_t iDigit=0; iDigit < cluster->GetNCells(); iDigit++)
1543 //Get from the absid the supermodule, tower and eta/phi numbers
1544 geom->GetCellIndex(cluster->GetCellAbsId(iDigit),iSupMod,iTower,iIphi,iIeta);
1545 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi,iIeta, iphi,ieta);
1547 //Get the cell energy, if recalibration is on, apply factors
1548 fraction = cluster->GetCellAmplitudeFraction(iDigit);
1549 if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1551 if (!fCellsRecalibrated)
1553 if(IsRecalibrationOn())
1555 recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1559 eCell = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1561 if(cluster->E() > 0 && eCell > 0)
1563 w = GetCellWeight(eCell,cluster->E());
1565 etai=(Double_t)ieta;
1566 phii=(Double_t)iphi;
1573 sEta += w * etai * etai ;
1574 etaMean += w * etai ;
1575 sPhi += w * phii * phii ;
1576 phiMean += w * phii ;
1577 sEtaPhi += 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;
1617 disp += w *((etai-etaMean)*(etai-etaMean)+(phii-phiMean)*(phii-phiMean));
1618 dEta += w * (etai-etaMean)*(etai-etaMean) ;
1619 dPhi += w * (phii-phiMean)*(phii-phiMean) ;
1623 AliError(Form("Wrong energy %f and/or amplitude %f\n", eCell, cluster->E()));
1626 //Normalize to the weigth and set shower shape parameters
1627 if (wtot > 0 && nstat > 1)
1636 sEta -= etaMean * etaMean ;
1637 sPhi -= phiMean * phiMean ;
1638 sEtaPhi -= etaMean * phiMean ;
1640 l0 = (0.5 * (sEta + sPhi) + TMath::Sqrt( 0.25 * (sEta - sPhi) * (sEta - sPhi) + sEtaPhi * sEtaPhi ));
1641 l1 = (0.5 * (sEta + sPhi) - TMath::Sqrt( 0.25 * (sEta - sPhi) * (sEta - sPhi) + sEtaPhi * sEtaPhi ));
1647 dEta = 0. ; dPhi = 0. ; disp = 0. ;
1648 sEta = 0. ; sPhi = 0. ; sEtaPhi = 0. ;
1653 //____________________________________________________________________________________________
1654 void AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(const AliEMCALGeometry * geom,
1655 AliVCaloCells* cells,
1656 AliVCluster * cluster)
1658 // Calculates new center of gravity in the local EMCAL-module coordinates
1659 // and tranfers into global ALICE coordinates
1660 // Calculates Dispersion and main axis and puts them into the cluster
1662 Float_t l0 = 0., l1 = 0.;
1663 Float_t disp = 0., dEta = 0., dPhi = 0.;
1664 Float_t sEta = 0., sPhi = 0., sEtaPhi = 0.;
1666 AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(geom,cells,cluster,l0,l1,disp,
1667 dEta, dPhi, sEta, sPhi, sEtaPhi);
1669 cluster->SetM02(l0);
1670 cluster->SetM20(l1);
1671 if(disp > 0. ) cluster->SetDispersion(TMath::Sqrt(disp)) ;
1675 //____________________________________________________________________________
1676 void AliEMCALRecoUtils::FindMatches(AliVEvent *event,
1677 TObjArray * clusterArr,
1678 const AliEMCALGeometry *geom)
1680 //This function should be called before the cluster loop
1681 //Before call this function, please recalculate the cluster positions
1682 //Given the input event, loop over all the tracks, select the closest cluster as matched with fCutR
1683 //Store matched cluster indexes and residuals
1685 fMatchedTrackIndex ->Reset();
1686 fMatchedClusterIndex->Reset();
1687 fResidualPhi->Reset();
1688 fResidualEta->Reset();
1690 fMatchedTrackIndex ->Set(1000);
1691 fMatchedClusterIndex->Set(1000);
1692 fResidualPhi->Set(1000);
1693 fResidualEta->Set(1000);
1695 AliESDEvent* esdevent = dynamic_cast<AliESDEvent*> (event);
1696 AliAODEvent* aodevent = dynamic_cast<AliAODEvent*> (event);
1698 // init the magnetic field if not already on
1699 if(!TGeoGlobalMagField::Instance()->GetField())
1701 AliInfo("Init the magnetic field\n");
1704 esdevent->InitMagneticField();
1708 Double_t curSol = 30000*aodevent->GetMagneticField()/5.00668;
1709 Double_t curDip = 6000 *aodevent->GetMuonMagFieldScale();
1710 AliMagF *field = AliMagF::CreateFieldMap(curSol,curDip);
1711 TGeoGlobalMagField::Instance()->SetField(field);
1715 AliInfo("Mag Field not initialized, null esd/aod evetn pointers");
1720 TObjArray *clusterArray = 0x0;
1723 clusterArray = new TObjArray(event->GetNumberOfCaloClusters());
1724 for(Int_t icl=0; icl<event->GetNumberOfCaloClusters(); icl++)
1726 AliVCluster *cluster = (AliVCluster*) event->GetCaloCluster(icl);
1727 if(geom && !IsGoodCluster(cluster,geom,(AliVCaloCells*)event->GetEMCALCells())) continue;
1728 clusterArray->AddAt(cluster,icl);
1734 for (Int_t i=0; i<21;i++) cv[i]=0;
1735 for(Int_t itr=0; itr<event->GetNumberOfTracks(); itr++)
1737 AliExternalTrackParam *trackParam = 0;
1739 //If the input event is ESD, the starting point for extrapolation is TPCOut, if available, or TPCInner
1740 AliESDtrack *esdTrack = 0;
1741 AliAODTrack *aodTrack = 0;
1744 esdTrack = esdevent->GetTrack(itr);
1745 if(!esdTrack) continue;
1746 if(!IsAccepted(esdTrack)) continue;
1747 if(esdTrack->Pt()<fCutMinTrackPt) continue;
1748 Double_t phi = esdTrack->Phi()*TMath::RadToDeg();
1749 if(TMath::Abs(esdTrack->Eta())>0.8 || phi <= 20 || phi >= 240 ) continue;
1750 trackParam = const_cast<AliExternalTrackParam*>(esdTrack->GetInnerParam());
1753 //If the input event is AOD, the starting point for extrapolation is at vertex
1754 //AOD tracks are selected according to its filterbit.
1757 aodTrack = aodevent->GetTrack(itr);
1758 if(!aodTrack) continue;
1759 if(!aodTrack->TestFilterMask(fAODFilterMask)) continue; //Select AOD tracks that fulfill GetStandardITSTPCTrackCuts2010()
1760 if(aodTrack->Pt()<fCutMinTrackPt) continue;
1761 Double_t phi = aodTrack->Phi()*TMath::RadToDeg();
1762 if(TMath::Abs(aodTrack->Eta())>0.8 || phi <= 20 || phi >= 240 ) continue;
1763 Double_t pos[3],mom[3];
1764 aodTrack->GetXYZ(pos);
1765 aodTrack->GetPxPyPz(mom);
1766 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()));
1767 trackParam= new AliExternalTrackParam(pos,mom,cv,aodTrack->Charge());
1770 //Return if the input data is not "AOD" or "ESD"
1773 printf("Wrong input data type! Should be \"AOD\" or \"ESD\"\n");
1776 clusterArray->Clear();
1777 delete clusterArray;
1782 if(!trackParam) continue;
1784 //Extrapolate the track to EMCal surface
1785 AliExternalTrackParam emcalParam(*trackParam);
1787 if(!ExtrapolateTrackToEMCalSurface(&emcalParam, 430., fMass, fStepSurface, eta, phi))
1789 if(aodevent && trackParam) delete trackParam;
1795 // esdTrack->SetOuterParam(&emcalParam,AliExternalTrackParam::kMultSec);
1798 if(TMath::Abs(eta)>0.75 || (phi) < 70*TMath::DegToRad() || (phi) > 190*TMath::DegToRad())
1800 if(aodevent && trackParam) delete trackParam;
1805 //Find matched clusters
1807 Float_t dEta = -999, dPhi = -999;
1810 index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArray, dEta, dPhi);
1814 index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArr, dEta, dPhi);
1819 fMatchedTrackIndex ->AddAt(itr,matched);
1820 fMatchedClusterIndex ->AddAt(index,matched);
1821 fResidualEta ->AddAt(dEta,matched);
1822 fResidualPhi ->AddAt(dPhi,matched);
1825 if(aodevent && trackParam) delete trackParam;
1830 clusterArray->Clear();
1831 delete clusterArray;
1834 AliDebug(2,Form("Number of matched pairs = %d !\n",matched));
1836 fMatchedTrackIndex ->Set(matched);
1837 fMatchedClusterIndex ->Set(matched);
1838 fResidualPhi ->Set(matched);
1839 fResidualEta ->Set(matched);
1842 //________________________________________________________________________________
1843 Int_t AliEMCALRecoUtils::FindMatchedClusterInEvent(const AliESDtrack *track,
1844 const AliVEvent *event,
1845 const AliEMCALGeometry *geom,
1846 Float_t &dEta, Float_t &dPhi)
1849 // This function returns the index of matched cluster to input track
1850 // Returns -1 if no match is found
1852 Double_t phiV = track->Phi()*TMath::RadToDeg();
1853 if(TMath::Abs(track->Eta())>0.8 || phiV <= 20 || phiV >= 240 ) return index;
1854 AliExternalTrackParam *trackParam = const_cast<AliExternalTrackParam*>(track->GetInnerParam());
1855 if(!trackParam) return index;
1856 AliExternalTrackParam emcalParam(*trackParam);
1858 if(!ExtrapolateTrackToEMCalSurface(&emcalParam, 430., fMass, fStepSurface, eta, phi)) return index;
1859 if(TMath::Abs(eta)>0.75 || (phi) < 70*TMath::DegToRad() || (phi) > 190*TMath::DegToRad()) return index;
1861 TObjArray *clusterArr = new TObjArray(event->GetNumberOfCaloClusters());
1863 for(Int_t icl=0; icl<event->GetNumberOfCaloClusters(); icl++)
1865 AliVCluster *cluster = (AliVCluster*) event->GetCaloCluster(icl);
1866 if(geom && !IsGoodCluster(cluster,geom,(AliVCaloCells*)event->GetEMCALCells())) continue;
1867 clusterArr->AddAt(cluster,icl);
1870 index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArr, dEta, dPhi);
1871 clusterArr->Clear();
1877 //_______________________________________________________________________________________________
1878 Int_t AliEMCALRecoUtils::FindMatchedClusterInClusterArr(const AliExternalTrackParam *emcalParam,
1879 AliExternalTrackParam *trkParam,
1880 const TObjArray * clusterArr,
1881 Float_t &dEta, Float_t &dPhi)
1883 // Find matched cluster in array
1885 dEta=-999, dPhi=-999;
1886 Float_t dRMax = fCutR, dEtaMax=fCutEta, dPhiMax=fCutPhi;
1888 Float_t tmpEta=-999, tmpPhi=-999;
1890 Double_t exPos[3] = {0.,0.,0.};
1891 if(!emcalParam->GetXYZ(exPos)) return index;
1893 Float_t clsPos[3] = {0.,0.,0.};
1894 for(Int_t icl=0; icl<clusterArr->GetEntriesFast(); icl++)
1896 AliVCluster *cluster = dynamic_cast<AliVCluster*> (clusterArr->At(icl)) ;
1897 if(!cluster || !cluster->IsEMCAL()) continue;
1898 cluster->GetPosition(clsPos);
1899 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));
1900 if(dR > fClusterWindow) continue;
1902 AliExternalTrackParam trkPamTmp (*trkParam);//Retrieve the starting point every time before the extrapolation
1903 if(!ExtrapolateTrackToCluster(&trkPamTmp, cluster, fMass, fStepCluster, tmpEta, tmpPhi)) continue;
1906 Float_t tmpR=TMath::Sqrt(tmpEta*tmpEta + tmpPhi*tmpPhi);
1915 else if(fCutEtaPhiSeparate)
1917 if(TMath::Abs(tmpEta)<TMath::Abs(dEtaMax) && TMath::Abs(tmpPhi)<TMath::Abs(dPhiMax))
1926 printf("Error: please specify your cut criteria\n");
1927 printf("To cut on sqrt(dEta^2+dPhi^2), use: SwitchOnCutEtaPhiSum()\n");
1928 printf("To cut on dEta and dPhi separately, use: SwitchOnCutEtaPhiSeparate()\n");
1939 //------------------------------------------------------------------------------------
1940 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToEMCalSurface(AliExternalTrackParam *trkParam,
1941 const Double_t emcalR,
1942 const Double_t mass,
1943 const Double_t step,
1947 //Extrapolate track to EMCAL surface
1949 eta = -999, phi = -999;
1950 if(!trkParam) return kFALSE;
1951 if(!AliTrackerBase::PropagateTrackToBxByBz(trkParam, emcalR, mass, step, kTRUE, 0.8, -1)) return kFALSE;
1952 Double_t trkPos[3] = {0.,0.,0.};
1953 if(!trkParam->GetXYZ(trkPos)) return kFALSE;
1954 TVector3 trkPosVec(trkPos[0],trkPos[1],trkPos[2]);
1955 eta = trkPosVec.Eta();
1956 phi = trkPosVec.Phi();
1958 phi += 2*TMath::Pi();
1963 //-----------------------------------------------------------------------------------
1964 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToPosition(AliExternalTrackParam *trkParam,
1965 const Float_t *clsPos,
1972 //Return the residual by extrapolating a track param to a global position
1976 if(!trkParam) return kFALSE;
1977 Double_t trkPos[3] = {0.,0.,0.};
1978 TVector3 vec(clsPos[0],clsPos[1],clsPos[2]);
1979 Double_t alpha = ((int)(vec.Phi()*TMath::RadToDeg()/20)+0.5)*20*TMath::DegToRad();
1980 vec.RotateZ(-alpha); //Rotate the cluster to the local extrapolation coordinate system
1981 if(!AliTrackerBase::PropagateTrackToBxByBz(trkParam, vec.X(), mass, step,kTRUE, 0.8, -1)) return kFALSE;
1982 if(!trkParam->GetXYZ(trkPos)) return kFALSE; //Get the extrapolated global position
1984 TVector3 clsPosVec(clsPos[0],clsPos[1],clsPos[2]);
1985 TVector3 trkPosVec(trkPos[0],trkPos[1],trkPos[2]);
1987 // track cluster matching
1988 tmpPhi = clsPosVec.DeltaPhi(trkPosVec); // tmpPhi is between -pi and pi
1989 tmpEta = clsPosVec.Eta()-trkPosVec.Eta();
1994 //----------------------------------------------------------------------------------
1995 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToCluster(AliExternalTrackParam *trkParam,
1996 const AliVCluster *cluster,
1997 const Double_t mass,
1998 const Double_t step,
2003 //Return the residual by extrapolating a track param to a cluster
2007 if(!cluster || !trkParam) return kFALSE;
2009 Float_t clsPos[3] = {0.,0.,0.};
2010 cluster->GetPosition(clsPos);
2012 return ExtrapolateTrackToPosition(trkParam, clsPos, mass, step, tmpEta, tmpPhi);
2015 //---------------------------------------------------------------------------------
2016 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToCluster(AliExternalTrackParam *trkParam,
2017 const AliVCluster *cluster,
2022 //Return the residual by extrapolating a track param to a clusterfStepCluster
2025 return ExtrapolateTrackToCluster(trkParam, cluster, fMass, fStepCluster, tmpEta, tmpPhi);
2028 //_______________________________________________________________________
2029 void AliEMCALRecoUtils::GetMatchedResiduals(const Int_t clsIndex,
2030 Float_t &dEta, Float_t &dPhi)
2032 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2033 //Get the residuals dEta and dPhi for this cluster to the closest track
2034 //Works with ESDs and AODs
2036 if( FindMatchedPosForCluster(clsIndex) >= 999 )
2038 AliDebug(2,"No matched tracks found!\n");
2043 dEta = fResidualEta->At(FindMatchedPosForCluster(clsIndex));
2044 dPhi = fResidualPhi->At(FindMatchedPosForCluster(clsIndex));
2047 //______________________________________________________________________________________________
2048 void AliEMCALRecoUtils::GetMatchedClusterResiduals(Int_t trkIndex, Float_t &dEta, Float_t &dPhi)
2050 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2051 //Get the residuals dEta and dPhi for this track to the closest cluster
2052 //Works with ESDs and AODs
2054 if( FindMatchedPosForTrack(trkIndex) >= 999 )
2056 AliDebug(2,"No matched cluster found!\n");
2061 dEta = fResidualEta->At(FindMatchedPosForTrack(trkIndex));
2062 dPhi = fResidualPhi->At(FindMatchedPosForTrack(trkIndex));
2065 //__________________________________________________________
2066 Int_t AliEMCALRecoUtils::GetMatchedTrackIndex(Int_t clsIndex)
2068 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2069 //Get the index of matched track to this cluster
2070 //Works with ESDs and AODs
2072 if(IsClusterMatched(clsIndex))
2073 return fMatchedTrackIndex->At(FindMatchedPosForCluster(clsIndex));
2078 //__________________________________________________________
2079 Int_t AliEMCALRecoUtils::GetMatchedClusterIndex(Int_t trkIndex)
2081 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2082 //Get the index of matched cluster to this track
2083 //Works with ESDs and AODs
2085 if(IsTrackMatched(trkIndex))
2086 return fMatchedClusterIndex->At(FindMatchedPosForTrack(trkIndex));
2091 //______________________________________________________________
2092 Bool_t AliEMCALRecoUtils::IsClusterMatched(Int_t clsIndex) const
2094 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2095 //Returns if the cluster has a match
2096 if(FindMatchedPosForCluster(clsIndex) < 999)
2102 //____________________________________________________________
2103 Bool_t AliEMCALRecoUtils::IsTrackMatched(Int_t trkIndex) const
2105 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2106 //Returns if the track has a match
2107 if(FindMatchedPosForTrack(trkIndex) < 999)
2113 //______________________________________________________________________
2114 UInt_t AliEMCALRecoUtils::FindMatchedPosForCluster(Int_t clsIndex) const
2116 //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2117 //Returns the position of the match in the fMatchedClusterIndex array
2118 Float_t tmpR = fCutR;
2121 for(Int_t i=0; i<fMatchedClusterIndex->GetSize(); i++)
2123 if(fMatchedClusterIndex->At(i)==clsIndex)
2125 Float_t r = TMath::Sqrt(fResidualEta->At(i)*fResidualEta->At(i) + fResidualPhi->At(i)*fResidualPhi->At(i));
2130 AliDebug(3,Form("Matched cluster index: index: %d, dEta: %2.4f, dPhi: %2.4f.\n",
2131 fMatchedClusterIndex->At(i),fResidualEta->At(i),fResidualPhi->At(i)));
2138 //____________________________________________________________________
2139 UInt_t AliEMCALRecoUtils::FindMatchedPosForTrack(Int_t trkIndex) const
2141 //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2142 //Returns the position of the match in the fMatchedTrackIndex array
2143 Float_t tmpR = fCutR;
2146 for(Int_t i=0; i<fMatchedTrackIndex->GetSize(); i++)
2148 if(fMatchedTrackIndex->At(i)==trkIndex)
2150 Float_t r = TMath::Sqrt(fResidualEta->At(i)*fResidualEta->At(i) + fResidualPhi->At(i)*fResidualPhi->At(i));
2155 AliDebug(3,Form("Matched track index: index: %d, dEta: %2.4f, dPhi: %2.4f.\n",
2156 fMatchedTrackIndex->At(i),fResidualEta->At(i),fResidualPhi->At(i)));
2163 //__________________________________________________________________________
2164 Bool_t AliEMCALRecoUtils::IsGoodCluster(AliVCluster *cluster,
2165 const AliEMCALGeometry *geom,
2166 AliVCaloCells* cells,const Int_t bc)
2168 // check if the cluster survives some quality cut
2171 Bool_t isGood=kTRUE;
2173 if(!cluster || !cluster->IsEMCAL()) return kFALSE;
2175 if(ClusterContainsBadChannel(geom,cluster->GetCellsAbsId(),cluster->GetNCells())) return kFALSE;
2177 if(!CheckCellFiducialRegion(geom,cluster,cells)) return kFALSE;
2179 if(IsExoticCluster(cluster, cells,bc)) return kFALSE;
2184 //__________________________________________________________
2185 Bool_t AliEMCALRecoUtils::IsAccepted(AliESDtrack *esdTrack)
2187 // Given a esd track, return whether the track survive all the cuts
2189 // The different quality parameter are first
2190 // retrieved from the track. then it is found out what cuts the
2191 // track did not survive and finally the cuts are imposed.
2193 UInt_t status = esdTrack->GetStatus();
2195 Int_t nClustersITS = esdTrack->GetITSclusters(0);
2196 Int_t nClustersTPC = esdTrack->GetTPCclusters(0);
2198 Float_t chi2PerClusterITS = -1;
2199 Float_t chi2PerClusterTPC = -1;
2200 if (nClustersITS!=0)
2201 chi2PerClusterITS = esdTrack->GetITSchi2()/Float_t(nClustersITS);
2202 if (nClustersTPC!=0)
2203 chi2PerClusterTPC = esdTrack->GetTPCchi2()/Float_t(nClustersTPC);
2207 if(fTrackCutsType==kGlobalCut)
2209 Float_t maxDCAToVertexXYPtDep = 0.0182 + 0.0350/TMath::Power(esdTrack->Pt(),1.01); //This expression comes from AliESDtrackCuts::GetStandardITSTPCTrackCuts2010()
2210 //AliDebug(3,Form("Track pT = %f, DCAtoVertexXY = %f",esdTrack->Pt(),MaxDCAToVertexXYPtDep));
2211 SetMaxDCAToVertexXY(maxDCAToVertexXYPtDep); //Set pT dependent DCA cut to vertex in x-y plane
2217 esdTrack->GetImpactParameters(b,bCov);
2218 if (bCov[0]<=0 || bCov[2]<=0)
2220 AliDebug(1, "Estimated b resolution lower or equal zero!");
2221 bCov[0]=0; bCov[2]=0;
2224 Float_t dcaToVertexXY = b[0];
2225 Float_t dcaToVertexZ = b[1];
2226 Float_t dcaToVertex = -1;
2228 if (fCutDCAToVertex2D)
2229 dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY/fCutMaxDCAToVertexXY/fCutMaxDCAToVertexXY + dcaToVertexZ*dcaToVertexZ/fCutMaxDCAToVertexZ/fCutMaxDCAToVertexZ);
2231 dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY + dcaToVertexZ*dcaToVertexZ);
2235 Bool_t cuts[kNCuts];
2236 for (Int_t i=0; i<kNCuts; i++) cuts[i]=kFALSE;
2238 // track quality cuts
2239 if (fCutRequireTPCRefit && (status&AliESDtrack::kTPCrefit)==0)
2241 if (fCutRequireITSRefit && (status&AliESDtrack::kITSrefit)==0)
2243 if (nClustersTPC<fCutMinNClusterTPC)
2245 if (nClustersITS<fCutMinNClusterITS)
2247 if (chi2PerClusterTPC>fCutMaxChi2PerClusterTPC)
2249 if (chi2PerClusterITS>fCutMaxChi2PerClusterITS)
2251 if (!fCutAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0)
2253 if (fCutDCAToVertex2D && dcaToVertex > 1)
2255 if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexXY) > fCutMaxDCAToVertexXY)
2257 if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexZ) > fCutMaxDCAToVertexZ)
2260 if(fTrackCutsType==kGlobalCut)
2262 //Require at least one SPD point + anything else in ITS
2263 if( (esdTrack->HasPointOnITSLayer(0) || esdTrack->HasPointOnITSLayer(1)) == kFALSE)
2268 for (Int_t i=0; i<kNCuts; i++)
2269 if (cuts[i]) { cut = kTRUE ; }
2278 //_____________________________________
2279 void AliEMCALRecoUtils::InitTrackCuts()
2281 //Intilize the track cut criteria
2282 //By default these cuts are set according to AliESDtrackCuts::GetStandardTPCOnlyTrackCuts()
2283 //Also you can customize the cuts using the setters
2285 switch (fTrackCutsType)
2289 AliInfo(Form("Track cuts for matching: GetStandardTPCOnlyTrackCuts()"));
2291 SetMinNClustersTPC(70);
2292 SetMaxChi2PerClusterTPC(4);
2293 SetAcceptKinkDaughters(kFALSE);
2294 SetRequireTPCRefit(kFALSE);
2297 SetRequireITSRefit(kFALSE);
2298 SetMaxDCAToVertexZ(3.2);
2299 SetMaxDCAToVertexXY(2.4);
2300 SetDCAToVertex2D(kTRUE);
2307 AliInfo(Form("Track cuts for matching: GetStandardITSTPCTrackCuts2010(kTURE)"));
2309 SetMinNClustersTPC(70);
2310 SetMaxChi2PerClusterTPC(4);
2311 SetAcceptKinkDaughters(kFALSE);
2312 SetRequireTPCRefit(kTRUE);
2315 SetRequireITSRefit(kTRUE);
2316 SetMaxDCAToVertexZ(2);
2317 SetMaxDCAToVertexXY();
2318 SetDCAToVertex2D(kFALSE);
2325 AliInfo(Form("Track cuts for matching: Loose cut w/o DCA cut"));
2326 SetMinNClustersTPC(50);
2327 SetAcceptKinkDaughters(kTRUE);
2335 //________________________________________________________________________
2336 void AliEMCALRecoUtils::SetClusterMatchedToTrack(const AliVEvent *event)
2338 // Checks if tracks are matched to EMC clusters and set the matched EMCAL cluster index to ESD track.
2340 Int_t nTracks = event->GetNumberOfTracks();
2341 for (Int_t iTrack = 0; iTrack < nTracks; ++iTrack)
2343 AliVTrack* track = dynamic_cast<AliVTrack*>(event->GetTrack(iTrack));
2346 AliWarning(Form("Could not receive track %d", iTrack));
2350 Int_t matchClusIndex = GetMatchedClusterIndex(iTrack);
2351 track->SetEMCALcluster(matchClusIndex); //sets -1 if track not matched within residual
2352 /*the following can be done better if AliVTrack::SetStatus will be there. Patch pending with Andreas/Peter*/
2353 AliESDtrack* esdtrack = dynamic_cast<AliESDtrack*>(track);
2355 if(matchClusIndex != -1)
2356 esdtrack->SetStatus(AliESDtrack::kEMCALmatch);
2358 esdtrack->ResetStatus(AliESDtrack::kEMCALmatch);
2360 AliAODTrack* aodtrack = dynamic_cast<AliAODTrack*>(track);
2361 if(matchClusIndex != -1)
2362 aodtrack->SetStatus(AliESDtrack::kEMCALmatch);
2364 aodtrack->ResetStatus(AliESDtrack::kEMCALmatch);
2368 AliDebug(2,"Track matched to closest cluster");
2371 //_________________________________________________________________________
2372 void AliEMCALRecoUtils::SetTracksMatchedToCluster(const AliVEvent *event)
2374 // Checks if EMC clusters are matched to ESD track.
2375 // Adds track indexes of all the tracks matched to a cluster withing residuals in ESDCalocluster.
2377 for (Int_t iClus=0; iClus < event->GetNumberOfCaloClusters(); ++iClus)
2379 AliVCluster *cluster = event->GetCaloCluster(iClus);
2380 if (!cluster->IsEMCAL())
2383 Int_t nTracks = event->GetNumberOfTracks();
2384 TArrayI arrayTrackMatched(nTracks);
2386 // Get the closest track matched to the cluster
2388 Int_t matchTrackIndex = GetMatchedTrackIndex(iClus);
2389 if (matchTrackIndex != -1)
2391 arrayTrackMatched[nMatched] = matchTrackIndex;
2395 // Get all other tracks matched to the cluster
2396 for(Int_t iTrk=0; iTrk<nTracks; ++iTrk)
2398 AliVTrack* track = dynamic_cast<AliVTrack*>(event->GetTrack(iTrk));
2399 if(iTrk == matchTrackIndex) continue;
2400 if(track->GetEMCALcluster() == iClus)
2402 arrayTrackMatched[nMatched] = iTrk;
2407 //printf("Tender::SetTracksMatchedToCluster - cluster E %f, N matches %d, first match %d\n",cluster->E(),nMatched,arrayTrackMatched[0]);
2409 arrayTrackMatched.Set(nMatched);
2410 AliESDCaloCluster *esdcluster = dynamic_cast<AliESDCaloCluster*>(cluster);
2412 esdcluster->AddTracksMatched(arrayTrackMatched);
2413 else if (nMatched>0) {
2414 AliAODCaloCluster *aodcluster = dynamic_cast<AliAODCaloCluster*>(cluster);
2416 aodcluster->AddTrackMatched(event->GetTrack(arrayTrackMatched.At(0)));
2419 Float_t eta= -999, phi = -999;
2420 if (matchTrackIndex != -1)
2421 GetMatchedResiduals(iClus, eta, phi);
2422 cluster->SetTrackDistance(phi, eta);
2425 AliDebug(2,"Cluster matched to tracks");
2428 //___________________________________________________
2429 void AliEMCALRecoUtils::Print(const Option_t *) const
2433 printf("AliEMCALRecoUtils Settings: \n");
2434 printf("Misalignment shifts\n");
2435 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,
2436 fMisalTransShift[i*3],fMisalTransShift[i*3+1],fMisalTransShift[i*3+2],
2437 fMisalRotShift[i*3], fMisalRotShift[i*3+1], fMisalRotShift[i*3+2] );
2438 printf("Non linearity function %d, parameters:\n", fNonLinearityFunction);
2439 for(Int_t i=0; i<6; i++) printf("param[%d]=%f\n",i, fNonLinearityParams[i]);
2441 printf("Position Recalculation option %d, Particle Type %d, fW0 %2.2f, Recalibrate Data %d \n",fPosAlgo,fParticleType,fW0, fRecalibration);
2443 printf("Matching criteria: ");
2446 printf("sqrt(dEta^2+dPhi^2)<%4.3f\n",fCutR);
2448 else if(fCutEtaPhiSeparate)
2450 printf("dEta<%4.3f, dPhi<%4.3f\n",fCutEta,fCutPhi);
2455 printf("please specify your cut criteria\n");
2456 printf("To cut on sqrt(dEta^2+dPhi^2), use: SwitchOnCutEtaPhiSum()\n");
2457 printf("To cut on dEta and dPhi separately, use: SwitchOnCutEtaPhiSeparate()\n");
2460 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);
2461 printf("Cluster selection window: dR < %2.0f\n",fClusterWindow);
2463 printf("Track cuts: \n");
2464 printf("Minimum track pT: %1.2f\n",fCutMinTrackPt);
2465 printf("AOD track selection mask: %d\n",fAODFilterMask);
2466 printf("TPCRefit = %d, ITSRefit = %d\n",fCutRequireTPCRefit,fCutRequireITSRefit);
2467 printf("AcceptKinks = %d\n",fCutAcceptKinkDaughters);
2468 printf("MinNCulsterTPC = %d, MinNClusterITS = %d\n",fCutMinNClusterTPC,fCutMinNClusterITS);
2469 printf("MaxChi2TPC = %2.2f, MaxChi2ITS = %2.2f\n",fCutMaxChi2PerClusterTPC,fCutMaxChi2PerClusterITS);
2470 printf("DCSToVertex2D = %d, MaxDCAToVertexXY = %2.2f, MaxDCAToVertexZ = %2.2f\n",fCutDCAToVertex2D,fCutMaxDCAToVertexXY,fCutMaxDCAToVertexZ);