]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALRecoUtils.cxx
df9784044082f151f05e49415a14c88f796c7976
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALRecoUtils.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /* $Id: AliEMCALRecoUtils.cxx 33808 2009-07-15 09:48:08Z gconesab $ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //
20 // Class AliEMCALRecoUtils
21 // Some utilities to recalculate the cluster position or energy linearity
22 //
23 //
24 // Author:  Gustavo Conesa (LPSC- Grenoble) 
25 //          Track matching part: Rongrong Ma (Yale)
26
27 ///////////////////////////////////////////////////////////////////////////////
28 // --- standard c ---
29
30 // standard C++ includes
31 //#include <Riostream.h>
32
33 // ROOT includes
34 #include <TGeoManager.h>
35 #include <TGeoMatrix.h>
36 #include <TGeoBBox.h>
37 #include <TH2F.h>
38 #include <TArrayI.h>
39 #include <TArrayF.h>
40 #include <TObjArray.h>
41
42 // STEER includes
43 #include "AliVCluster.h"
44 #include "AliVCaloCells.h"
45 #include "AliLog.h"
46 #include "AliPID.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"
54
55 // EMCAL includes
56 #include "AliEMCALRecoUtils.h"
57 #include "AliEMCALGeometry.h"
58 #include "AliTrackerBase.h"
59 #include "AliEMCALPIDUtils.h"
60
61
62 ClassImp(AliEMCALRecoUtils)
63   
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)
85 {
86 //
87   // Constructor.
88   // Initialize all constant values which have to be used
89   // during Reco algorithm execution
90   //
91   
92   // Init parameters
93   InitParameters();
94   
95   //Track matching
96   fMatchedTrackIndex     = new TArrayI();
97   fMatchedClusterIndex   = new TArrayI();
98   fResidualPhi           = new TArrayF();
99   fResidualEta           = new TArrayF();
100   fPIDUtils              = new AliEMCALPIDUtils();
101
102   InitTrackCuts();
103 }
104
105 //______________________________________________________________________
106 AliEMCALRecoUtils::AliEMCALRecoUtils(const AliEMCALRecoUtils & reco) 
107 : TNamed(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)
136 {
137   //Copy ctor
138   
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]  ; }
143
144 }
145
146
147 //______________________________________________________________________
148 AliEMCALRecoUtils & AliEMCALRecoUtils::operator = (const AliEMCALRecoUtils & reco) 
149 {
150   //Assignment operator
151   
152   if(this == &reco)return *this;
153   ((TNamed *)this)->operator=(reco);
154
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]  ; }   
159   
160   fParticleType              = reco.fParticleType;
161   fPosAlgo                   = reco.fPosAlgo; 
162   fW0                        = reco.fW0;
163   
164   fNonLinearityFunction      = reco.fNonLinearityFunction;
165   fNonLinearThreshold        = reco.fNonLinearThreshold;
166   fSmearClusterEnergy        = reco.fSmearClusterEnergy;
167
168   fCellsRecalibrated         = reco.fCellsRecalibrated;
169   fRecalibration             = reco.fRecalibration;
170   fEMCALRecalibrationFactors = reco.fEMCALRecalibrationFactors;
171
172   fTimeRecalibration             = reco.fTimeRecalibration;
173   fEMCALTimeRecalibrationFactors = reco.fEMCALTimeRecalibrationFactors;
174
175   fUseRunCorrectionFactors   = reco.fUseRunCorrectionFactors;
176   
177   fRemoveBadChannels         = reco.fRemoveBadChannels;
178   fRecalDistToBadChannels    = reco.fRecalDistToBadChannels;
179   fEMCALBadChannelMap        = reco.fEMCALBadChannelMap;
180   
181   fNCellsFromEMCALBorder     = reco.fNCellsFromEMCALBorder;
182   fNoEMCALBorderAtEta0       = reco.fNoEMCALBorderAtEta0;
183   
184   fRejectExoticCluster       = reco.fRejectExoticCluster;           
185   fRejectExoticCells         = reco.fRejectExoticCells; 
186   fExoticCellFraction        = reco.fExoticCellFraction;
187   fExoticCellDiffTime        = reco.fExoticCellDiffTime;              
188   fExoticCellMinAmplitude    = reco.fExoticCellMinAmplitude;
189   
190   fPIDUtils                  = reco.fPIDUtils;
191
192   fAODFilterMask             = reco.fAODFilterMask;
193   
194   fCutEtaPhiSum              = reco.fCutEtaPhiSum;
195   fCutEtaPhiSeparate         = reco.fCutEtaPhiSeparate;
196   fCutR                      = reco.fCutR;
197   fCutEta                    = reco.fCutEta;
198   fCutPhi                    = reco.fCutPhi;
199   fClusterWindow             = reco.fClusterWindow;
200   fMass                      = reco.fMass;
201   fStepSurface               = reco.fStepSurface;
202   fStepCluster               = reco.fStepCluster;
203
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;
216   
217   if(reco.fResidualEta)
218   {
219     // assign or copy construct
220     if(fResidualEta)
221     { 
222       *fResidualEta = *reco.fResidualEta;
223     }
224     else 
225     {
226       fResidualEta = new TArrayF(*reco.fResidualEta);
227     }
228   }
229   else
230   {
231     if(fResidualEta)delete fResidualEta;
232     fResidualEta = 0;
233   }
234   
235   if(reco.fResidualPhi)
236   {
237     // assign or copy construct
238     if(fResidualPhi)
239     { 
240       *fResidualPhi = *reco.fResidualPhi;
241     }
242     else 
243     {
244       fResidualPhi = new TArrayF(*reco.fResidualPhi);
245     }
246   }
247   else
248   {
249     if(fResidualPhi)delete fResidualPhi;
250     fResidualPhi = 0;
251   }
252   
253   if(reco.fMatchedTrackIndex)
254   {
255     // assign or copy construct
256     if(fMatchedTrackIndex)
257     { 
258       *fMatchedTrackIndex = *reco.fMatchedTrackIndex;
259     }
260     else 
261     { 
262       fMatchedTrackIndex = new TArrayI(*reco.fMatchedTrackIndex);
263     }
264   }
265   else
266   {
267     if(fMatchedTrackIndex)delete fMatchedTrackIndex;
268     fMatchedTrackIndex = 0;
269   }  
270   
271   if(reco.fMatchedClusterIndex)
272   {
273     // assign or copy construct
274     if(fMatchedClusterIndex)
275     { 
276       *fMatchedClusterIndex = *reco.fMatchedClusterIndex;
277     }
278     else 
279     {
280       fMatchedClusterIndex = new TArrayI(*reco.fMatchedClusterIndex);
281     }
282   }
283   else
284   {
285     if(fMatchedClusterIndex)delete fMatchedClusterIndex;
286     fMatchedClusterIndex = 0;
287   }
288    
289   return *this;
290 }
291
292
293 //_____________________________________
294 AliEMCALRecoUtils::~AliEMCALRecoUtils()
295 {
296   //Destructor.
297   
298   if(fEMCALRecalibrationFactors) 
299   { 
300     fEMCALRecalibrationFactors->Clear();
301     delete fEMCALRecalibrationFactors;
302   }  
303   
304   if(fEMCALTimeRecalibrationFactors) 
305   { 
306     fEMCALTimeRecalibrationFactors->Clear();
307     delete fEMCALTimeRecalibrationFactors;
308   }  
309   
310   if(fEMCALBadChannelMap) 
311   { 
312     fEMCALBadChannelMap->Clear();
313     delete fEMCALBadChannelMap;
314   }
315  
316   delete fMatchedTrackIndex   ; 
317   delete fMatchedClusterIndex ; 
318   delete fResidualEta         ; 
319   delete fResidualPhi         ; 
320   delete fPIDUtils            ;
321
322   InitTrackCuts();
323 }
324
325 //_______________________________________________________________________________
326 Bool_t AliEMCALRecoUtils::AcceptCalibrateCell(const Int_t absID, const Int_t bc,
327                                               Float_t  & amp,    Double_t & time, 
328                                               AliVCaloCells* cells) 
329 {
330   // Reject cell if criteria not passed and calibrate it
331   
332   AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
333   
334   if(absID < 0 || absID >= 24*48*geom->GetNumberOfSuperModules()) return kFALSE;
335   
336   Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1; 
337   
338   if(!geom->GetCellIndex(absID,imod,iTower,iIphi,iIeta)) 
339   {
340     // cell absID does not exist
341     amp=0; time = 1.e9;
342     return kFALSE; 
343   }
344   
345   geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);  
346
347   // Do not include bad channels found in analysis,
348   if( IsBadChannelsRemovalSwitchedOn() && GetEMCALChannelStatus(imod, ieta, iphi)) 
349   {
350     return kFALSE;
351   }
352   
353   //Recalibrate energy
354   amp  = cells->GetCellAmplitude(absID);
355   if(!fCellsRecalibrated && IsRecalibrationOn())
356     amp *= GetEMCALChannelRecalibrationFactor(imod,ieta,iphi);
357   
358   
359   // Recalibrate time
360   time = cells->GetCellTime(absID);
361   
362   RecalibrateCellTime(absID,bc,time);
363   
364   return kTRUE;
365 }
366
367 //_____________________________________________________________________________
368 Bool_t AliEMCALRecoUtils::CheckCellFiducialRegion(const AliEMCALGeometry* geom, 
369                                                   const AliVCluster* cluster, 
370                                                   AliVCaloCells* cells) 
371 {
372   // Given the list of AbsId of the cluster, get the maximum cell and 
373   // check if there are fNCellsFromBorder from the calorimeter border
374   
375   if(!cluster)
376   {
377     AliInfo("Cluster pointer null!");
378     return kFALSE;
379   }
380   
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;
383   
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);
387   
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));
390   
391   if(absIdMax==-1) return kFALSE;
392   
393   //Check if the cell is close to the borders:
394   Bool_t okrow = kFALSE;
395   Bool_t okcol = kFALSE;
396   
397   if(iSM < 0 || iphi < 0 || ieta < 0 ) 
398   {
399     AliFatal(Form("Negative value for super module: %d, or cell ieta: %d, or cell iphi: %d, check EMCAL geometry name\n",
400                   iSM,ieta,iphi));
401   }
402   
403   //Check rows/phi
404   if(iSM < 10)
405   {
406     if(iphi >= fNCellsFromEMCALBorder && iphi < 24-fNCellsFromEMCALBorder) okrow =kTRUE; 
407   }
408   else if (iSM >=10 && ( ( geom->GetEMCGeometry()->GetGeoName()).Contains("12SMV1"))) 
409   {
410     if(iphi >= fNCellsFromEMCALBorder && iphi < 8-fNCellsFromEMCALBorder) okrow =kTRUE; //1/3 sm case
411   }
412   else 
413   {
414     if(iphi >= fNCellsFromEMCALBorder && iphi < 12-fNCellsFromEMCALBorder) okrow =kTRUE; // half SM case
415   }
416   
417   //Check columns/eta
418   if(!fNoEMCALBorderAtEta0)
419   {
420     if(ieta  > fNCellsFromEMCALBorder && ieta < 48-fNCellsFromEMCALBorder) okcol =kTRUE; 
421   }
422   else
423   {
424     if(iSM%2==0)
425     {
426       if(ieta >= fNCellsFromEMCALBorder)     okcol = kTRUE;  
427     }
428     else 
429     {
430       if(ieta <  48-fNCellsFromEMCALBorder)  okcol = kTRUE;  
431     }
432   }//eta 0 not checked
433   
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));
436   
437   if (okcol && okrow) 
438   {
439     //printf("Accept\n");
440     return kTRUE;
441   }
442   else  
443   {
444     //printf("Reject\n");
445     AliDebug(2,Form("Reject cluster in border, max cell : ieta %d, iphi %d, SM %d\n",ieta, iphi, iSM));
446     return kFALSE;
447   }
448   
449 }  
450
451
452 //_______________________________________________________________________________
453 Bool_t AliEMCALRecoUtils::ClusterContainsBadChannel(const AliEMCALGeometry* geom, 
454                                                     const UShort_t* cellList, 
455                                                     const Int_t nCells)
456 {
457   // Check that in the cluster cells, there is no bad channel of those stored 
458   // in fEMCALBadChannelMap or fPHOSBadChannelMap
459   
460   if(!fRemoveBadChannels)  return kFALSE;
461   if(!fEMCALBadChannelMap) return kFALSE;
462   
463   Int_t icol = -1;
464   Int_t irow = -1;
465   Int_t imod = -1;
466   for(Int_t iCell = 0; iCell<nCells; iCell++)
467   {
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))
474     {
475       AliDebug(2,Form("Cluster with bad channel: SM %d, col %d, row %d\n",imod, icol, irow));
476       return kTRUE;
477     }
478     
479   }// cell cluster loop
480   
481   return kFALSE;
482 }
483
484
485 //___________________________________________________________________________
486 Float_t AliEMCALRecoUtils::GetECross(const Int_t absID, const Double_t tcell,
487                                      AliVCaloCells* cells, const Int_t bc)
488 {
489   //Calculate the energy in the cross around the energy given cell
490   
491   AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
492   
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);  
496   
497   //Get close cells index, energy and time, not in corners
498   
499   Int_t absID1 = -1;
500   Int_t absID2 = -1;
501   
502   if( iphi < AliEMCALGeoParams::fgkEMCALRows-1) absID1 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi+1, ieta);
503   if( iphi > 0 )                                absID2 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi-1, ieta);
504   
505   // In case of cell in eta = 0 border, depending on SM shift the cross cell index
506   
507   Int_t absID3 = -1;
508   Int_t absID4 = -1;
509   
510   if     ( ieta == AliEMCALGeoParams::fgkEMCALCols-1 && !(imod%2) )
511   {
512     absID3 = geom-> GetAbsCellIdFromCellIndexes(imod+1, iphi, 0);
513     absID4 = geom-> GetAbsCellIdFromCellIndexes(imod,   iphi, ieta-1); 
514   }
515   else if( ieta == 0 && imod%2 )
516   {
517     absID3 = geom-> GetAbsCellIdFromCellIndexes(imod,   iphi, ieta+1);
518     absID4 = geom-> GetAbsCellIdFromCellIndexes(imod-1, iphi, AliEMCALGeoParams::fgkEMCALCols-1); 
519   }
520   else
521   {
522     if( ieta < AliEMCALGeoParams::fgkEMCALCols-1 ) 
523       absID3 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta+1);
524     if( ieta > 0 )                                 
525       absID4 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta-1); 
526   }
527   
528   //printf("IMOD %d, AbsId %d, a %d, b %d, c %d e %d \n",imod,absID,absID1,absID2,absID3,absID4);
529   
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;
533   
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); 
538   
539   /*
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);
550    */
551   
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 ;
556   
557   return ecell1+ecell2+ecell3+ecell4;
558   
559 }
560
561 //_____________________________________________________________________________________________
562 Bool_t AliEMCALRecoUtils::IsExoticCell(const Int_t absID, AliVCaloCells* cells, const Int_t bc)
563 {
564   // Look to cell neighbourhood and reject if it seems exotic
565   // Do before recalibrating the cells
566
567   if(!fRejectExoticCells) return kFALSE;
568   
569   Float_t  ecell  = 0;
570   Double_t tcell  = 0;
571   Bool_t   accept = AcceptCalibrateCell(absID, bc, ecell ,tcell ,cells); 
572   
573   if(!accept) return kTRUE; // reject this cell
574   
575   if(ecell < fExoticCellMinAmplitude) return kFALSE; // do not reject low energy cells
576
577   Float_t eCross = GetECross(absID,tcell,cells,bc);
578   
579   if(1-eCross/ecell > fExoticCellFraction) 
580   {
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));
583     return kTRUE;
584   }
585
586   return kFALSE;
587 }
588
589 //___________________________________________________________________
590 Bool_t AliEMCALRecoUtils::IsExoticCluster(const AliVCluster *cluster, 
591                                           AliVCaloCells *cells, 
592                                           const Int_t bc) 
593 {
594   // Check if the cluster highest energy tower is exotic
595   
596   if(!cluster)
597   {
598     AliInfo("Cluster pointer null!");
599     return kFALSE;
600   }
601   
602   if(!fRejectExoticCluster) return kFALSE;
603   
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);
609   
610   return IsExoticCell(absId,cells,bc);
611   
612 }
613
614 //_______________________________________________________________________
615 Float_t AliEMCALRecoUtils::SmearClusterEnergy(const AliVCluster* cluster) 
616 {
617   //In case of MC analysis, smear energy to match resolution/calibration in real data
618   
619   if(!cluster)
620   {
621     AliInfo("Cluster pointer null!");
622     return 0;
623   }
624   
625   Float_t energy    = cluster->E() ;
626   Float_t rdmEnergy = energy ;
627   if(fSmearClusterEnergy)
628   {
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));
633   }
634   
635   return rdmEnergy;
636 }
637
638 //____________________________________________________________________________
639 Float_t AliEMCALRecoUtils::CorrectClusterEnergyLinearity(AliVCluster* cluster)
640 {
641   // Correct cluster energy from non linearity functions
642   
643   if(!cluster)
644   {
645     AliInfo("Cluster pointer null!");
646     return 0;
647   }
648   
649   Float_t energy = cluster->E();
650
651   switch (fNonLinearityFunction) 
652   {
653       
654     case kPi0MC:
655     {
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]))));
665       break;
666     }
667      
668     case kPi0MCv2:
669     {
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;      
674       
675       energy *= fNonLinearityParams[0]/TMath::Power(energy+fNonLinearityParams[1],fNonLinearityParams[2])+1;
676       break;
677     }
678       
679     case kPi0GammaGamma:
680     {
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
686       break;
687     }
688       
689     case kPi0GammaConversion:
690     {
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));
696       
697       break;
698     }
699       
700     case kBeamTest:
701     {
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]));
709       
710       break;
711     }
712       
713     case kBeamTestCorrected:
714     {
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]))));
724
725       break;
726     }
727       
728     case kNoCorrection:
729       AliDebug(2,"No correction on the energy\n");
730       break;
731       
732   }
733
734   return energy;
735 }
736
737 //__________________________________________________
738 void AliEMCALRecoUtils::InitNonLinearityParam()
739 {
740   //Initialising Non Linearity Parameters
741   
742   if(fNonLinearityFunction == kPi0MC) 
743   {
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;
749   }
750   
751   if(fNonLinearityFunction == kPi0GammaGamma) 
752   {
753     fNonLinearityParams[0] = 1.04;
754     fNonLinearityParams[1] = -0.1445;
755     fNonLinearityParams[2] = 1.046;
756   }  
757
758   if(fNonLinearityFunction == kPi0GammaConversion) 
759   {
760     fNonLinearityParams[0] = 0.139393;
761     fNonLinearityParams[1] = 0.0566186;
762     fNonLinearityParams[2] = 0.982133;
763   }  
764
765   if(fNonLinearityFunction == kBeamTest) 
766   {
767     if(fNonLinearThreshold == 30) 
768     {
769       fNonLinearityParams[0] = 1.007; 
770       fNonLinearityParams[1] = 0.894; 
771       fNonLinearityParams[2] = 0.246; 
772     }
773     if(fNonLinearThreshold == 45) 
774     {
775       fNonLinearityParams[0] = 1.003; 
776       fNonLinearityParams[1] = 0.719; 
777       fNonLinearityParams[2] = 0.334; 
778     }
779     if(fNonLinearThreshold == 75) 
780     {
781       fNonLinearityParams[0] = 1.002; 
782       fNonLinearityParams[1] = 0.797; 
783       fNonLinearityParams[2] = 0.358; 
784     }
785   }
786
787   if(fNonLinearityFunction == kBeamTestCorrected) 
788   {
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;
796   }
797 }
798
799 //_________________________________________________________
800 Float_t  AliEMCALRecoUtils::GetDepth(const Float_t energy, 
801                                      const Int_t iParticle, 
802                                      const Int_t iSM) const 
803 {
804   //Calculate shower depth for a given cluster energy and particle type
805
806   // parameters 
807   Float_t x0    = 1.31;
808   Float_t ecr   = 8;
809   Float_t depth = 0;
810   
811   switch ( iParticle )
812   {
813     case kPhoton:
814       depth = x0 * (TMath::Log(energy*1000/ ecr) + 0.5); //Multiply energy by 1000 to transform to MeV
815       break;
816       
817     case kElectron:
818       depth = x0 * (TMath::Log(energy*1000/ ecr) - 0.5); //Multiply energy by 1000 to transform to MeV
819       break;
820       
821     case kHadron:
822       // hadron 
823       // boxes anc. here
824       if(gGeoManager)
825       {
826         gGeoManager->cd("ALIC_1/XEN1_1");
827         TGeoNode        *geoXEn1    = gGeoManager->GetCurrentNode();
828         TGeoNodeMatrix  *geoSM      = dynamic_cast<TGeoNodeMatrix *>(geoXEn1->GetDaughter(iSM));
829         if(geoSM)
830         {
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");
836         }
837         else AliFatal("NULL  GEANT node matrix");
838       }
839       else
840       {//electron
841         depth = x0 * (TMath::Log(energy*1000 / ecr)  - 0.5); //Multiply energy by 1000 to transform to MeV
842       }
843         
844       break;
845       
846     default://photon
847       depth = x0 * (TMath::Log(energy*1000 / ecr) + 0.5); //Multiply energy by 1000 to transform to MeV
848   }  
849   
850   return depth;
851 }
852
853 //____________________________________________________________________
854 void AliEMCALRecoUtils::GetMaxEnergyCell(const AliEMCALGeometry *geom, 
855                                          AliVCaloCells* cells, 
856                                          const AliVCluster* clu, 
857                                          Int_t  & absId,  
858                                          Int_t  & iSupMod, 
859                                          Int_t  & ieta, 
860                                          Int_t  & iphi, 
861                                          Bool_t & shared)
862 {
863   //For a given CaloCluster gets the absId of the cell 
864   //with maximum energy deposit.
865   
866   Double_t eMax        = -1.;
867   Double_t eCell       = -1.;
868   Float_t  fraction    = 1.;
869   Float_t  recalFactor = 1.;
870   Int_t    cellAbsId   = -1 ;
871
872   Int_t iTower  = -1;
873   Int_t iIphi   = -1;
874   Int_t iIeta   = -1;
875   Int_t iSupMod0= -1;
876
877   if(!clu)
878   {
879     AliInfo("Cluster pointer null!");
880     absId=-1; iSupMod0=-1, ieta = -1; iphi = -1; shared = -1;
881     return;
882   }
883   
884   for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++) 
885   {
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);
892     if     (iDig==0) 
893     {
894       iSupMod0=iSupMod;
895     }
896     else if(iSupMod0!=iSupMod) 
897     {
898       shared = kTRUE;
899       //printf("AliEMCALRecoUtils::GetMaxEnergyCell() - SHARED CLUSTER\n");
900     }
901     if(!fCellsRecalibrated && IsRecalibrationOn()) 
902     {
903       recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
904     }
905     eCell  = cells->GetCellAmplitude(cellAbsId)*fraction*recalFactor;
906     //printf("b Cell %d, id, %d, amp %f, fraction %f\n",iDig,cellAbsId,eCell,fraction);
907     if(eCell > eMax)  
908     { 
909       eMax  = eCell; 
910       absId = cellAbsId;
911       //printf("\t new max: cell %d, e %f, ecell %f\n",maxId, eMax,eCell);
912     }
913   }// cell loop
914   
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");
922 }
923
924 //______________________________________
925 void AliEMCALRecoUtils::InitParameters()
926 {
927   // Initialize data members with default values
928   
929   fParticleType = kPhoton;
930   fPosAlgo      = kUnchanged;
931   fW0           = 4.5;
932   
933   fNonLinearityFunction = kNoCorrection;
934   fNonLinearThreshold   = 30;
935   
936   fExoticCellFraction     = 0.97;
937   fExoticCellDiffTime     = 1e6;
938   fExoticCellMinAmplitude = 0.5;
939   
940   fAODFilterMask = 32;
941   
942   fCutEtaPhiSum      = kTRUE;
943   fCutEtaPhiSeparate = kFALSE;
944   
945   fCutR   = 0.05; 
946   fCutEta = 0.025; 
947   fCutPhi = 0.05;
948   
949   fClusterWindow = 100;
950   fMass          = 0.139;
951   
952   fStepSurface   = 20.;                      
953   fStepCluster   = 5.;
954   fTrackCutsType = kLooseCut;
955   
956   fCutMinTrackPt     = 0;
957   fCutMinNClusterTPC = -1;
958   fCutMinNClusterITS = -1;
959   
960   fCutMaxChi2PerClusterTPC  = 1e10;
961   fCutMaxChi2PerClusterITS  = 1e10;
962   
963   fCutRequireTPCRefit     = kFALSE;
964   fCutRequireITSRefit     = kFALSE;
965   fCutAcceptKinkDaughters = kFALSE;
966   
967   fCutMaxDCAToVertexXY = 1e10;             
968   fCutMaxDCAToVertexZ  = 1e10;              
969   fCutDCAToVertex2D    = kFALSE;
970   
971   
972   //Misalignment matrices
973   for(Int_t i = 0; i < 15 ; i++) 
974   {
975     fMisalTransShift[i] = 0.; 
976     fMisalRotShift[i]   = 0.; 
977   }
978   
979   //Non linearity
980   for(Int_t i = 0; i < 7  ; i++) fNonLinearityParams[i] = 0.; 
981   
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;
990   
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;
995   
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
1001 }
1002
1003 //_____________________________________________________
1004 void AliEMCALRecoUtils::InitEMCALRecalibrationFactors()
1005 {
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);
1011   
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++) 
1018   {
1019     for (Int_t i = 0; i < 48; i++) 
1020     {
1021       for (Int_t j = 0; j < 24; j++) 
1022       {
1023         SetEMCALChannelRecalibrationFactor(sm,i,j,1.);
1024       }
1025     }
1026   }
1027   
1028   fEMCALRecalibrationFactors->SetOwner(kTRUE);
1029   fEMCALRecalibrationFactors->Compress();
1030   
1031   //In order to avoid rewriting the same histograms
1032   TH1::AddDirectory(oldStatus);    
1033 }
1034
1035 //_________________________________________________________
1036 void AliEMCALRecoUtils::InitEMCALTimeRecalibrationFactors()
1037 {
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);
1043   
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++) 
1051   {
1052     for (Int_t i = 0; i < 48*24*12; i++) 
1053       SetEMCALChannelTimeRecalibrationFactor(bc,i,0.);
1054   }
1055   
1056   fEMCALTimeRecalibrationFactors->SetOwner(kTRUE);
1057   fEMCALTimeRecalibrationFactors->Compress();
1058   
1059   //In order to avoid rewriting the same histograms
1060   TH1::AddDirectory(oldStatus);    
1061 }
1062
1063 //____________________________________________________
1064 void AliEMCALRecoUtils::InitEMCALBadChannelStatusMap()
1065 {
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);
1071   
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++) 
1075   {
1076     fEMCALBadChannelMap->Add(new TH2I(Form("EMCALBadChannelMap_Mod%d",i),Form("EMCALBadChannelMap_Mod%d",i), 48, 0, 48, 24, 0, 24));
1077   }
1078   
1079   fEMCALBadChannelMap->SetOwner(kTRUE);
1080   fEMCALBadChannelMap->Compress();
1081   
1082   //In order to avoid rewriting the same histograms
1083   TH1::AddDirectory(oldStatus);    
1084 }
1085
1086 //____________________________________________________________________________
1087 void AliEMCALRecoUtils::RecalibrateClusterEnergy(const AliEMCALGeometry* geom, 
1088                                                  AliVCluster * cluster, 
1089                                                  AliVCaloCells * cells, 
1090                                                  const Int_t bc)
1091 {
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();
1095   
1096   if(!cluster)
1097   {
1098     AliInfo("Cluster pointer null!");
1099     return;
1100   }  
1101   
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();
1106   
1107   //Initialize some used variables
1108   Float_t energy = 0;
1109   Int_t   absId  =-1;
1110   Int_t   icol   =-1, irow =-1, imod=1;
1111   Float_t factor = 1, frac = 0;
1112   Int_t   absIdMax = -1;
1113   Float_t emax     = 0;
1114   
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++)
1117   {
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
1121     
1122     if(!fCellsRecalibrated && IsRecalibrationOn()) 
1123     {
1124       // Energy  
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);
1130       
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)));
1133       
1134     } 
1135     
1136     energy += cells->GetCellAmplitude(absId)*factor*frac;
1137     
1138     if(emax < cells->GetCellAmplitude(absId)*factor*frac)
1139     {
1140       emax     = cells->GetCellAmplitude(absId)*factor*frac;
1141       absIdMax = absId;
1142     }
1143   }
1144   
1145   AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - Energy before %f, after %f \n",cluster->E(),energy));
1146
1147   cluster->SetE(energy);
1148
1149   // Recalculate time of cluster   
1150   Double_t timeorg = cluster->GetTOF();
1151   if(!fCellsRecalibrated && IsTimeRecalibrationOn())
1152   {
1153     Double_t time = timeorg;
1154     RecalibrateCellTime(absIdMax,bc,time);
1155     cluster->SetTOF(time);
1156   } 
1157
1158   AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - Time before %f, after %f \n",timeorg,cluster->GetTOF()));
1159
1160 }
1161
1162 //_____________________________________________________________
1163 void AliEMCALRecoUtils::RecalibrateCells(AliVCaloCells * cells,
1164                                          const Int_t bc)
1165 {
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();
1169
1170   if(!IsRecalibrationOn() && !IsTimeRecalibrationOn()) return;
1171   
1172   if(!cells)
1173   {
1174     AliInfo("Cells pointer null!");
1175     return;
1176   }  
1177   
1178   Short_t  absId  =-1;
1179   Bool_t   accept = kFALSE;
1180   Float_t  ecell  = 0;
1181   Double_t tcell  = 0;
1182   Double_t ecellin = 0;
1183   Double_t tcellin = 0;
1184   Short_t  mclabel = -1;
1185   Double_t efrac = 0;
1186   
1187   Int_t nEMcell  = cells->GetNumberOfCells() ;  
1188   for (Int_t iCell = 0; iCell < nEMcell; iCell++) 
1189   { 
1190     cells->GetCell( iCell, absId, ecellin, tcellin, mclabel, efrac );
1191     
1192     accept = AcceptCalibrateCell(absId, bc, ecell ,tcell ,cells); 
1193     if(!accept) 
1194     {
1195       ecell = 0;
1196       tcell = -1;
1197     }
1198     
1199     //Set new values
1200     cells->SetCell(iCell,absId,ecell, tcell, mclabel, efrac);
1201   }
1202
1203   fCellsRecalibrated = kTRUE;
1204 }
1205
1206 //_______________________________________________________________________________________________________
1207 void AliEMCALRecoUtils::RecalibrateCellTime(const Int_t absId, const Int_t bc, Double_t & celltime) const
1208 {
1209   // Recalibrate time of cell with absID  considering the recalibration map 
1210   // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1211   
1212   if(!fCellsRecalibrated && IsTimeRecalibrationOn() && bc >= 0)
1213   {
1214     celltime -= GetEMCALChannelTimeRecalibrationFactor(bc%4,absId)*1.e-9;    ;  
1215   }
1216 }
1217   
1218 //______________________________________________________________________________
1219 void AliEMCALRecoUtils::RecalculateClusterPosition(const AliEMCALGeometry *geom, 
1220                                                    AliVCaloCells* cells, 
1221                                                    AliVCluster* clu)
1222 {
1223   //For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1224   
1225   if(!clu)
1226   {
1227     AliInfo("Cluster pointer null!");
1228     return;
1229   }
1230     
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.");
1234 }  
1235
1236 //_____________________________________________________________________________________________
1237 void AliEMCALRecoUtils::RecalculateClusterPositionFromTowerGlobal(const AliEMCALGeometry *geom, 
1238                                                                   AliVCaloCells* cells, 
1239                                                                   AliVCluster* clu)
1240 {
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
1243   
1244   Double_t eCell       = 0.;
1245   Float_t  fraction    = 1.;
1246   Float_t  recalFactor = 1.;
1247   
1248   Int_t    absId   = -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;
1255
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) ;
1259   
1260   //printf("** Cluster energy %f, ncells %d, depth %f\n",clEnergy,clu->GetNCells(),depth);
1261   
1262   for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++) 
1263   {
1264     absId = clu->GetCellAbsId(iDig);
1265     fraction  = clu->GetCellAmplitudeFraction(iDig);
1266     if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1267     
1268     if (!fCellsRecalibrated)
1269     {
1270       geom->GetCellIndex(absId,iSM,iTower,iIphi,iIeta); 
1271       geom->GetCellPhiEtaIndexInSModule(iSM,iTower,iIphi, iIeta,iphi,ieta);      
1272       
1273       if(IsRecalibrationOn()) 
1274       {
1275         recalFactor = GetEMCALChannelRecalibrationFactor(iSM,ieta,iphi);
1276       }
1277     }
1278     
1279     eCell  = cells->GetCellAmplitude(absId)*fraction*recalFactor;
1280     
1281     weight = GetCellWeight(eCell,clEnergy);
1282     totalWeight += weight;
1283     
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]);
1288
1289     for(int i=0; i<3; i++ ) newPos[i] += (weight*pGlobal[i]);
1290   }// cell loop
1291   
1292   if(totalWeight>0)
1293   {
1294     for(int i=0; i<3; i++ )    newPos[i] /= totalWeight;
1295   }
1296     
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]);
1301   
1302   if(iSupModMax > 1) //sector 1
1303   {
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]);
1308   } else //sector 0
1309   {
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]);
1314   }
1315   //printf("NewPos : %2.3f,%2.3f,%2.3f\n",newPos[0],newPos[1],newPos[2]);
1316
1317   clu->SetPosition(newPos);
1318 }  
1319
1320 //____________________________________________________________________________________________
1321 void AliEMCALRecoUtils::RecalculateClusterPositionFromTowerIndex(const AliEMCALGeometry *geom, 
1322                                                                  AliVCaloCells* cells, 
1323                                                                  AliVCluster* clu)
1324 {
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
1327   
1328   Double_t eCell       = 1.;
1329   Float_t  fraction    = 1.;
1330   Float_t  recalFactor = 1.;
1331   
1332   Int_t absId   = -1;
1333   Int_t iTower  = -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;
1338
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) ;
1342
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;
1346   
1347   for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++) 
1348   {
1349     absId = clu->GetCellAbsId(iDig);
1350     fraction  = clu->GetCellAmplitudeFraction(iDig);
1351     if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1352
1353     if     (iDig==0)  startingSM = iSupMod;
1354     else if(iSupMod != startingSM) areInSameSM = kFALSE;
1355
1356     eCell  = cells->GetCellAmplitude(absId);
1357     
1358     geom->GetCellIndex(absId,iSupMod,iTower,iIphi,iIeta); 
1359     geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);    
1360     
1361     if (!fCellsRecalibrated)
1362     {
1363       if(IsRecalibrationOn()) 
1364       {
1365         recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1366       }
1367     }
1368     
1369     eCell  = cells->GetCellAmplitude(absId)*fraction*recalFactor;
1370     
1371     weight = GetCellWeight(eCell,clEnergy);
1372     if(weight < 0) weight = 0;
1373     totalWeight += weight;
1374     weightedCol += ieta*weight;
1375     weightedRow += iphi*weight;
1376     
1377     //printf("Max cell? cell %d, amplitude org %f, fraction %f, recalibration %f, amplitude new %f \n",cellAbsId, cells->GetCellAmplitude(cellAbsId), fraction, recalFactor, eCell) ;
1378   }// cell loop
1379     
1380   Float_t xyzNew[]={0.,0.,0.};
1381   if(areInSameSM == kTRUE) 
1382   {
1383     //printf("In Same SM\n");
1384     weightedCol = weightedCol/totalWeight;
1385     weightedRow = weightedRow/totalWeight;
1386     geom->RecalculateTowerPosition(weightedRow, weightedCol, iSupModMax, depth, fMisalTransShift, fMisalRotShift, xyzNew); 
1387   } 
1388   else 
1389   {
1390     //printf("In Different SM\n");
1391     geom->RecalculateTowerPosition(iphi,        ieta,        iSupModMax, depth, fMisalTransShift, fMisalRotShift, xyzNew); 
1392   }
1393   
1394   clu->SetPosition(xyzNew);
1395 }
1396
1397 //___________________________________________________________________________________________
1398 void AliEMCALRecoUtils::RecalculateClusterDistanceToBadChannel(const AliEMCALGeometry * geom, 
1399                                                                AliVCaloCells* cells, 
1400                                                                AliVCluster * cluster)
1401 {           
1402   //re-evaluate distance to bad channel with updated bad map
1403   
1404   if(!fRecalDistToBadChannels) return;
1405   
1406   if(!cluster)
1407   {
1408     AliInfo("Cluster pointer null!");
1409     return;
1410   }  
1411   
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);
1417
1418   Int_t dRrow, dRcol;  
1419   Float_t  minDist = 10000.;
1420   Float_t  dist    = 0.;
1421   
1422   //Loop on tower status map 
1423   for(Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++)
1424   {
1425     for(Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++)
1426     {
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);
1431       
1432       dRrow=TMath::Abs(irowM-irow);
1433       dRcol=TMath::Abs(icolM-icol);
1434       dist=TMath::Sqrt(dRrow*dRrow+dRcol*dRcol);
1435       if(dist < minDist)
1436       {
1437         //printf("MIN DISTANCE TO BAD %2.2f\n",dist);
1438         minDist = dist;
1439       }
1440     }
1441   }
1442   
1443   //In case the cluster is shared by 2 SuperModules, need to check the map of the second Super Module
1444   if (shared) 
1445   {
1446     TH2D* hMap2 = 0;
1447     Int_t iSupMod2 = -1;
1448     
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);
1453     
1454     //Loop on tower status map of second super module
1455     for(Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++)
1456     {
1457       for(Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++)
1458       {
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);
1464         
1465         if(iSupMod%2) 
1466         {
1467           dRcol=TMath::Abs(icol-(AliEMCALGeoParams::fgkEMCALCols+icolM));
1468         } else 
1469         {
1470           dRcol=TMath::Abs(AliEMCALGeoParams::fgkEMCALCols+icol-icolM);
1471         }                    
1472         
1473         dist=TMath::Sqrt(dRrow*dRrow+dRcol*dRcol);
1474         if(dist < minDist) minDist = dist;        
1475       }
1476     }
1477   }// shared cluster in 2 SuperModules
1478   
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);
1481 }
1482
1483 //__________________________________________________________________
1484 void AliEMCALRecoUtils::RecalculateClusterPID(AliVCluster * cluster)
1485 {           
1486   //re-evaluate identification parameters with bayesian
1487   
1488   if(!cluster)
1489   {
1490     AliInfo("Cluster pointer null!");
1491     return;
1492   }
1493   
1494   if ( cluster->GetM02() != 0)
1495     fPIDUtils->ComputePID(cluster->E(),cluster->GetM02());
1496   
1497   Float_t pidlist[AliPID::kSPECIESN+1];
1498   for(Int_t i = 0; i < AliPID::kSPECIESN+1; i++) pidlist[i] = fPIDUtils->GetPIDFinal(i);
1499         
1500   cluster->SetPID(pidlist);
1501 }
1502
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)
1510 {
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
1514   
1515   if(!cluster)
1516   {
1517     AliInfo("Cluster pointer null!");
1518     return;
1519   }
1520     
1521   Double_t eCell       = 0.;
1522   Float_t  fraction    = 1.;
1523   Float_t  recalFactor = 1.;
1524
1525   Int_t    iSupMod = -1;
1526   Int_t    iTower  = -1;
1527   Int_t    iIphi   = -1;
1528   Int_t    iIeta   = -1;
1529   Int_t    iphi    = -1;
1530   Int_t    ieta    = -1;
1531   Double_t etai    = -1.;
1532   Double_t phii    = -1.;
1533   
1534   Int_t    nstat   = 0 ;
1535   Float_t  wtot    = 0.;
1536   Double_t w       = 0.;
1537   Double_t etaMean = 0.;
1538   Double_t phiMean = 0.;
1539     
1540   //Loop on cells
1541   for(Int_t iDigit=0; iDigit < cluster->GetNCells(); iDigit++) 
1542   {
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);        
1546     
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
1550     
1551     if (!fCellsRecalibrated)
1552     {
1553       if(IsRecalibrationOn()) 
1554       {
1555         recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1556       }
1557     }
1558     
1559     eCell  = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1560     
1561     if(cluster->E() > 0 && eCell > 0)
1562     {
1563       w  = GetCellWeight(eCell,cluster->E());
1564       
1565       etai=(Double_t)ieta;
1566       phii=(Double_t)iphi;  
1567       
1568       if(w > 0.0) 
1569       {
1570         wtot += w ;
1571         nstat++;            
1572         //Shower shape
1573         sEta     += w * etai * etai ;
1574         etaMean  += w * etai ;
1575         sPhi     += w * phii * phii ;
1576         phiMean  += w * phii ; 
1577         sEtaPhi  += w * etai * phii ; 
1578       }
1579     } 
1580     else
1581       AliError(Form("Wrong energy %f and/or amplitude %f\n", eCell, cluster->E()));
1582   }//cell loop
1583   
1584   //Normalize to the weight  
1585   if (wtot > 0) 
1586   {
1587     etaMean /= wtot ;
1588     phiMean /= wtot ;
1589   }
1590   else
1591     AliError(Form("Wrong weight %f\n", wtot));
1592   
1593   //Calculate dispersion  
1594   for(Int_t iDigit=0; iDigit < cluster->GetNCells(); iDigit++) 
1595   {
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);
1599     
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()) 
1604     {
1605       recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1606     }
1607     eCell  = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1608     
1609     if(cluster->E() > 0 && eCell > 0)
1610     {
1611       w  = GetCellWeight(eCell,cluster->E());
1612       
1613       etai=(Double_t)ieta;
1614       phii=(Double_t)iphi;    
1615       if(w > 0.0) 
1616       { 
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) ; 
1620       }
1621     }
1622     else
1623       AliError(Form("Wrong energy %f and/or amplitude %f\n", eCell, cluster->E()));
1624   }// cell loop
1625   
1626   //Normalize to the weigth and set shower shape parameters
1627   if (wtot > 0 && nstat > 1) 
1628   {
1629     disp    /= wtot ;
1630     dEta    /= wtot ;
1631     dPhi    /= wtot ;
1632     sEta    /= wtot ;
1633     sPhi    /= wtot ;
1634     sEtaPhi /= wtot ;
1635     
1636     sEta    -= etaMean * etaMean ;
1637     sPhi    -= phiMean * phiMean ;
1638     sEtaPhi -= etaMean * phiMean ;
1639     
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 ));
1642   }
1643   else
1644   {
1645     l0   = 0. ;
1646     l1   = 0. ;
1647     dEta = 0. ; dPhi = 0. ; disp    = 0. ;
1648     sEta = 0. ; sPhi = 0. ; sEtaPhi = 0. ;
1649   }  
1650   
1651 }
1652
1653 //____________________________________________________________________________________________
1654 void AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(const AliEMCALGeometry * geom, 
1655                                                                 AliVCaloCells* cells, 
1656                                                                 AliVCluster * cluster)
1657 {
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
1661   
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.;
1665   
1666   AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(geom,cells,cluster,l0,l1,disp,
1667                                                              dEta, dPhi, sEta, sPhi, sEtaPhi);
1668   
1669   cluster->SetM02(l0);
1670   cluster->SetM20(l1);
1671   if(disp > 0. ) cluster->SetDispersion(TMath::Sqrt(disp)) ;
1672   
1673
1674
1675 //____________________________________________________________________________
1676 void AliEMCALRecoUtils::FindMatches(AliVEvent *event,
1677                                     TObjArray * clusterArr,  
1678                                     const AliEMCALGeometry *geom)
1679 {
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
1684   
1685   fMatchedTrackIndex  ->Reset();
1686   fMatchedClusterIndex->Reset();
1687   fResidualPhi->Reset();
1688   fResidualEta->Reset();
1689   
1690   fMatchedTrackIndex  ->Set(1000);
1691   fMatchedClusterIndex->Set(1000);
1692   fResidualPhi->Set(1000);
1693   fResidualEta->Set(1000);
1694   
1695   AliESDEvent* esdevent = dynamic_cast<AliESDEvent*> (event);
1696   AliAODEvent* aodevent = dynamic_cast<AliAODEvent*> (event);
1697   
1698   // init the magnetic field if not already on
1699   if(!TGeoGlobalMagField::Instance()->GetField())
1700   {
1701     AliInfo("Init the magnetic field\n");
1702     if     (esdevent) 
1703     {
1704       esdevent->InitMagneticField();
1705     }
1706     else if(aodevent)
1707     {
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);
1712     }
1713     else
1714     {
1715       AliInfo("Mag Field not initialized, null esd/aod evetn pointers");
1716     }
1717     
1718   } // Init mag field
1719   
1720   TObjArray *clusterArray = 0x0;
1721   if(!clusterArr)
1722     {
1723       clusterArray = new TObjArray(event->GetNumberOfCaloClusters());
1724       for(Int_t icl=0; icl<event->GetNumberOfCaloClusters(); icl++)
1725   {
1726     AliVCluster *cluster = (AliVCluster*) event->GetCaloCluster(icl);
1727     if(geom && !IsGoodCluster(cluster,geom,(AliVCaloCells*)event->GetEMCALCells())) continue;
1728     clusterArray->AddAt(cluster,icl);
1729   }
1730     }
1731   
1732   Int_t    matched=0;
1733   Double_t cv[21];
1734   for (Int_t i=0; i<21;i++) cv[i]=0;
1735   for(Int_t itr=0; itr<event->GetNumberOfTracks(); itr++)
1736   {
1737     AliExternalTrackParam *trackParam = 0;
1738     
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;
1742     if(esdevent)
1743     {
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());
1751     }
1752     
1753     //If the input event is AOD, the starting point for extrapolation is at vertex
1754     //AOD tracks are selected according to its filterbit.
1755     else if(aodevent)
1756     {
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());
1768     }
1769     
1770     //Return if the input data is not "AOD" or "ESD"
1771     else
1772     {
1773       printf("Wrong input data type! Should be \"AOD\" or \"ESD\"\n");
1774       if(clusterArray)
1775   {
1776     clusterArray->Clear();
1777     delete clusterArray;
1778   }
1779       return;
1780     }
1781     
1782     if(!trackParam) continue;
1783
1784     //Extrapolate the track to EMCal surface
1785     AliExternalTrackParam emcalParam(*trackParam);
1786     Float_t eta, phi;
1787     if(!ExtrapolateTrackToEMCalSurface(&emcalParam, 430., fMass, fStepSurface, eta, phi)) 
1788       {
1789   if(aodevent && trackParam) delete trackParam;
1790   continue;
1791       }
1792
1793 //    if(esdevent)
1794 //      {
1795 //  esdTrack->SetOuterParam(&emcalParam,AliExternalTrackParam::kMultSec);
1796 //      }
1797
1798     if(TMath::Abs(eta)>0.75 || (phi) < 70*TMath::DegToRad() || (phi) > 190*TMath::DegToRad())
1799       {
1800   if(aodevent && trackParam) delete trackParam;
1801   continue;
1802       }
1803
1804
1805     //Find matched clusters
1806     Int_t index = -1;
1807     Float_t dEta = -999, dPhi = -999;
1808     if(!clusterArr)
1809       {
1810   index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArray, dEta, dPhi);  
1811       }
1812     else
1813       {
1814   index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArr, dEta, dPhi);  
1815       }  
1816     
1817     if(index>-1)
1818     {
1819       fMatchedTrackIndex   ->AddAt(itr,matched);
1820       fMatchedClusterIndex ->AddAt(index,matched);
1821       fResidualEta         ->AddAt(dEta,matched);
1822       fResidualPhi         ->AddAt(dPhi,matched);
1823       matched++;
1824     }
1825     if(aodevent && trackParam) delete trackParam;
1826   }//track loop
1827
1828   if(clusterArray)
1829     {
1830       clusterArray->Clear();
1831       delete clusterArray;
1832     }
1833   
1834   AliDebug(2,Form("Number of matched pairs = %d !\n",matched));
1835   
1836   fMatchedTrackIndex   ->Set(matched);
1837   fMatchedClusterIndex ->Set(matched);
1838   fResidualPhi         ->Set(matched);
1839   fResidualEta         ->Set(matched);
1840 }
1841
1842 //________________________________________________________________________________
1843 Int_t AliEMCALRecoUtils::FindMatchedClusterInEvent(const AliESDtrack *track, 
1844                                                    const AliVEvent *event, 
1845                                                    const AliEMCALGeometry *geom, 
1846                                                    Float_t &dEta, Float_t &dPhi)
1847 {
1848   //
1849   // This function returns the index of matched cluster to input track
1850   // Returns -1 if no match is found
1851   Int_t index = -1;
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);
1857   Float_t eta, phi;
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;
1860
1861   TObjArray *clusterArr = new TObjArray(event->GetNumberOfCaloClusters());
1862
1863   for(Int_t icl=0; icl<event->GetNumberOfCaloClusters(); icl++)
1864   {
1865     AliVCluster *cluster = (AliVCluster*) event->GetCaloCluster(icl);
1866     if(geom && !IsGoodCluster(cluster,geom,(AliVCaloCells*)event->GetEMCALCells())) continue;
1867     clusterArr->AddAt(cluster,icl);
1868   }
1869
1870   index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArr, dEta, dPhi);  
1871   clusterArr->Clear();
1872   delete clusterArr;
1873   
1874   return index;
1875 }
1876
1877 //_______________________________________________________________________________________________
1878 Int_t  AliEMCALRecoUtils::FindMatchedClusterInClusterArr(const AliExternalTrackParam *emcalParam, 
1879                                                          AliExternalTrackParam *trkParam, 
1880                                                          const TObjArray * clusterArr, 
1881                                                          Float_t &dEta, Float_t &dPhi)
1882 {
1883   // Find matched cluster in array
1884   
1885   dEta=-999, dPhi=-999;
1886   Float_t dRMax = fCutR, dEtaMax=fCutEta, dPhiMax=fCutPhi;
1887   Int_t index = -1;
1888   Float_t tmpEta=-999, tmpPhi=-999;
1889
1890   Double_t exPos[3] = {0.,0.,0.};
1891   if(!emcalParam->GetXYZ(exPos)) return index;
1892
1893   Float_t clsPos[3] = {0.,0.,0.};
1894   for(Int_t icl=0; icl<clusterArr->GetEntriesFast(); icl++)
1895     {
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;
1901
1902       AliExternalTrackParam trkPamTmp (*trkParam);//Retrieve the starting point every time before the extrapolation
1903       if(!ExtrapolateTrackToCluster(&trkPamTmp, cluster, fMass, fStepCluster, tmpEta, tmpPhi)) continue;
1904       if(fCutEtaPhiSum)
1905         {
1906           Float_t tmpR=TMath::Sqrt(tmpEta*tmpEta + tmpPhi*tmpPhi);
1907           if(tmpR<dRMax)
1908       {
1909         dRMax=tmpR;
1910         dEtaMax=tmpEta;
1911         dPhiMax=tmpPhi;
1912         index=icl;
1913       }
1914         }
1915       else if(fCutEtaPhiSeparate)
1916         {
1917           if(TMath::Abs(tmpEta)<TMath::Abs(dEtaMax) && TMath::Abs(tmpPhi)<TMath::Abs(dPhiMax))
1918       {
1919         dEtaMax = tmpEta;
1920         dPhiMax = tmpPhi;
1921         index=icl;
1922       }
1923         }
1924       else
1925         {
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");
1929           return index;
1930         }
1931     }
1932
1933   dEta=dEtaMax;
1934   dPhi=dPhiMax;
1935
1936   return index;
1937 }
1938
1939 //------------------------------------------------------------------------------------
1940 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToEMCalSurface(AliExternalTrackParam *trkParam, 
1941                                                          const Double_t emcalR,
1942                                                          const Double_t mass, 
1943                                                          const Double_t step, 
1944                                                          Float_t &eta, 
1945                                                          Float_t &phi)
1946 {
1947   //Extrapolate track to EMCAL surface
1948   
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();
1957   if(phi<0)
1958     phi += 2*TMath::Pi();
1959
1960   return kTRUE;
1961 }
1962
1963 //-----------------------------------------------------------------------------------
1964 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToPosition(AliExternalTrackParam *trkParam, 
1965                                                      const Float_t *clsPos, 
1966                                                      Double_t mass, 
1967                                                      Double_t step, 
1968                                                      Float_t &tmpEta, 
1969                                                      Float_t &tmpPhi)
1970 {
1971   //
1972   //Return the residual by extrapolating a track param to a global position
1973   //
1974   tmpEta = -999;
1975   tmpPhi = -999;
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
1983
1984   TVector3 clsPosVec(clsPos[0],clsPos[1],clsPos[2]);
1985   TVector3 trkPosVec(trkPos[0],trkPos[1],trkPos[2]);
1986
1987   // track cluster matching
1988   tmpPhi = clsPosVec.DeltaPhi(trkPosVec);    // tmpPhi is between -pi and pi
1989   tmpEta = clsPosVec.Eta()-trkPosVec.Eta();
1990
1991   return kTRUE;
1992 }
1993
1994 //----------------------------------------------------------------------------------
1995 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToCluster(AliExternalTrackParam *trkParam, 
1996                                                     const AliVCluster *cluster, 
1997                                                     const Double_t mass, 
1998                                                     const Double_t step, 
1999                                                     Float_t &tmpEta, 
2000                                                     Float_t &tmpPhi)
2001 {
2002   //
2003   //Return the residual by extrapolating a track param to a cluster
2004   //
2005   tmpEta = -999;
2006   tmpPhi = -999;
2007   if(!cluster || !trkParam) return kFALSE;
2008
2009   Float_t clsPos[3] = {0.,0.,0.};
2010   cluster->GetPosition(clsPos);
2011
2012   return ExtrapolateTrackToPosition(trkParam, clsPos, mass, step, tmpEta, tmpPhi);
2013 }
2014
2015 //---------------------------------------------------------------------------------
2016 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToCluster(AliExternalTrackParam *trkParam, 
2017                                                     const AliVCluster *cluster, 
2018                                                     Float_t &tmpEta, 
2019                                                     Float_t &tmpPhi)
2020 {
2021   //
2022   //Return the residual by extrapolating a track param to a clusterfStepCluster
2023   //
2024
2025   return ExtrapolateTrackToCluster(trkParam, cluster, fMass, fStepCluster, tmpEta, tmpPhi);
2026 }
2027
2028 //_______________________________________________________________________
2029 void AliEMCALRecoUtils::GetMatchedResiduals(const Int_t clsIndex, 
2030                                             Float_t &dEta, Float_t &dPhi)
2031 {
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
2035
2036   if( FindMatchedPosForCluster(clsIndex) >= 999 )
2037   {
2038     AliDebug(2,"No matched tracks found!\n");
2039     dEta=999.;
2040     dPhi=999.;
2041     return;
2042   }
2043   dEta = fResidualEta->At(FindMatchedPosForCluster(clsIndex));
2044   dPhi = fResidualPhi->At(FindMatchedPosForCluster(clsIndex));
2045 }
2046
2047 //______________________________________________________________________________________________
2048 void AliEMCALRecoUtils::GetMatchedClusterResiduals(Int_t trkIndex, Float_t &dEta, Float_t &dPhi)
2049 {
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
2053
2054   if( FindMatchedPosForTrack(trkIndex) >= 999 )
2055   {
2056     AliDebug(2,"No matched cluster found!\n");
2057     dEta=999.;
2058     dPhi=999.;
2059     return;
2060   }
2061   dEta = fResidualEta->At(FindMatchedPosForTrack(trkIndex));
2062   dPhi = fResidualPhi->At(FindMatchedPosForTrack(trkIndex));
2063 }
2064
2065 //__________________________________________________________
2066 Int_t AliEMCALRecoUtils::GetMatchedTrackIndex(Int_t clsIndex)
2067 {
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
2071   
2072   if(IsClusterMatched(clsIndex))
2073     return fMatchedTrackIndex->At(FindMatchedPosForCluster(clsIndex));
2074   else 
2075     return -1; 
2076 }
2077
2078 //__________________________________________________________
2079 Int_t AliEMCALRecoUtils::GetMatchedClusterIndex(Int_t trkIndex)
2080 {
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
2084   
2085   if(IsTrackMatched(trkIndex))
2086     return fMatchedClusterIndex->At(FindMatchedPosForTrack(trkIndex));
2087   else 
2088     return -1; 
2089 }
2090
2091 //______________________________________________________________
2092 Bool_t AliEMCALRecoUtils::IsClusterMatched(Int_t clsIndex) const
2093 {
2094   //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2095   //Returns if the cluster has a match
2096   if(FindMatchedPosForCluster(clsIndex) < 999) 
2097     return kTRUE;
2098   else
2099     return kFALSE;
2100 }
2101
2102 //____________________________________________________________
2103 Bool_t AliEMCALRecoUtils::IsTrackMatched(Int_t trkIndex) const 
2104 {
2105   //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2106   //Returns if the track has a match
2107   if(FindMatchedPosForTrack(trkIndex) < 999) 
2108     return kTRUE;
2109   else
2110     return kFALSE;
2111 }
2112
2113 //______________________________________________________________________
2114 UInt_t AliEMCALRecoUtils::FindMatchedPosForCluster(Int_t clsIndex) const
2115 {
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;
2119   UInt_t pos = 999;
2120   
2121   for(Int_t i=0; i<fMatchedClusterIndex->GetSize(); i++) 
2122   {
2123     if(fMatchedClusterIndex->At(i)==clsIndex) 
2124     {
2125       Float_t r = TMath::Sqrt(fResidualEta->At(i)*fResidualEta->At(i) + fResidualPhi->At(i)*fResidualPhi->At(i));
2126       if(r<tmpR) 
2127       {
2128         pos=i;
2129         tmpR=r;
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)));
2132       }
2133     }
2134   }
2135   return pos;
2136 }
2137
2138 //____________________________________________________________________
2139 UInt_t AliEMCALRecoUtils::FindMatchedPosForTrack(Int_t trkIndex) const
2140 {
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;
2144   UInt_t pos = 999;
2145   
2146   for(Int_t i=0; i<fMatchedTrackIndex->GetSize(); i++) 
2147   {
2148     if(fMatchedTrackIndex->At(i)==trkIndex) 
2149     {
2150       Float_t r = TMath::Sqrt(fResidualEta->At(i)*fResidualEta->At(i) + fResidualPhi->At(i)*fResidualPhi->At(i));
2151       if(r<tmpR) 
2152       {
2153         pos=i;
2154         tmpR=r;
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)));
2157       }
2158     }
2159   }
2160   return pos;
2161 }
2162
2163 //__________________________________________________________________________
2164 Bool_t AliEMCALRecoUtils::IsGoodCluster(AliVCluster *cluster, 
2165                                         const AliEMCALGeometry *geom, 
2166                                         AliVCaloCells* cells,const Int_t bc)
2167 {
2168   // check if the cluster survives some quality cut
2169   //
2170   //
2171   Bool_t isGood=kTRUE;
2172
2173   if(!cluster || !cluster->IsEMCAL())              return kFALSE;
2174   
2175   if(ClusterContainsBadChannel(geom,cluster->GetCellsAbsId(),cluster->GetNCells())) return kFALSE;
2176   
2177   if(!CheckCellFiducialRegion(geom,cluster,cells)) return kFALSE;
2178   
2179   if(IsExoticCluster(cluster, cells,bc))           return kFALSE;
2180
2181   return isGood;
2182 }
2183
2184 //__________________________________________________________
2185 Bool_t AliEMCALRecoUtils::IsAccepted(AliESDtrack *esdTrack)
2186 {
2187   // Given a esd track, return whether the track survive all the cuts
2188
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.
2192
2193   UInt_t status = esdTrack->GetStatus();
2194
2195   Int_t nClustersITS = esdTrack->GetITSclusters(0);
2196   Int_t nClustersTPC = esdTrack->GetTPCclusters(0);
2197
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);
2204
2205
2206   //DCA cuts
2207   if(fTrackCutsType==kGlobalCut)
2208     {
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
2212     }
2213
2214
2215   Float_t b[2];
2216   Float_t bCov[3];
2217   esdTrack->GetImpactParameters(b,bCov);
2218   if (bCov[0]<=0 || bCov[2]<=0) 
2219   {
2220     AliDebug(1, "Estimated b resolution lower or equal zero!");
2221     bCov[0]=0; bCov[2]=0;
2222   }
2223
2224   Float_t dcaToVertexXY = b[0];
2225   Float_t dcaToVertexZ = b[1];
2226   Float_t dcaToVertex = -1;
2227
2228   if (fCutDCAToVertex2D)
2229     dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY/fCutMaxDCAToVertexXY/fCutMaxDCAToVertexXY + dcaToVertexZ*dcaToVertexZ/fCutMaxDCAToVertexZ/fCutMaxDCAToVertexZ);
2230   else
2231     dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY + dcaToVertexZ*dcaToVertexZ);
2232     
2233   // cut the track?
2234   
2235   Bool_t cuts[kNCuts];
2236   for (Int_t i=0; i<kNCuts; i++) cuts[i]=kFALSE;
2237   
2238   // track quality cuts
2239   if (fCutRequireTPCRefit && (status&AliESDtrack::kTPCrefit)==0)
2240     cuts[0]=kTRUE;
2241   if (fCutRequireITSRefit && (status&AliESDtrack::kITSrefit)==0)
2242     cuts[1]=kTRUE;
2243   if (nClustersTPC<fCutMinNClusterTPC)
2244     cuts[2]=kTRUE;
2245   if (nClustersITS<fCutMinNClusterITS) 
2246     cuts[3]=kTRUE;
2247   if (chi2PerClusterTPC>fCutMaxChi2PerClusterTPC) 
2248     cuts[4]=kTRUE; 
2249   if (chi2PerClusterITS>fCutMaxChi2PerClusterITS) 
2250     cuts[5]=kTRUE;  
2251   if (!fCutAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0)
2252     cuts[6]=kTRUE;
2253   if (fCutDCAToVertex2D && dcaToVertex > 1)
2254     cuts[7] = kTRUE;
2255   if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexXY) > fCutMaxDCAToVertexXY)
2256     cuts[8] = kTRUE;
2257   if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexZ) > fCutMaxDCAToVertexZ)
2258     cuts[9] = kTRUE;
2259
2260   if(fTrackCutsType==kGlobalCut)
2261     {
2262       //Require at least one SPD point + anything else in ITS
2263       if( (esdTrack->HasPointOnITSLayer(0) || esdTrack->HasPointOnITSLayer(1)) == kFALSE)
2264   cuts[10] = kTRUE;
2265     }
2266
2267   Bool_t cut=kFALSE;
2268   for (Int_t i=0; i<kNCuts; i++)
2269     if (cuts[i]) { cut = kTRUE ; }
2270
2271     // cut the track
2272   if (cut) 
2273     return kFALSE;
2274   else 
2275     return kTRUE;
2276 }
2277
2278 //_____________________________________
2279 void AliEMCALRecoUtils::InitTrackCuts()
2280 {
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
2284   
2285   switch (fTrackCutsType)
2286   {
2287     case kTPCOnlyCut:
2288     {
2289       AliInfo(Form("Track cuts for matching: GetStandardTPCOnlyTrackCuts()"));
2290       //TPC
2291       SetMinNClustersTPC(70);
2292       SetMaxChi2PerClusterTPC(4);
2293       SetAcceptKinkDaughters(kFALSE);
2294       SetRequireTPCRefit(kFALSE);
2295       
2296       //ITS
2297       SetRequireITSRefit(kFALSE);
2298       SetMaxDCAToVertexZ(3.2);
2299       SetMaxDCAToVertexXY(2.4);
2300       SetDCAToVertex2D(kTRUE);
2301       
2302       break;
2303     }
2304       
2305     case kGlobalCut:
2306     {
2307       AliInfo(Form("Track cuts for matching: GetStandardITSTPCTrackCuts2010(kTURE)"));
2308       //TPC
2309       SetMinNClustersTPC(70);
2310       SetMaxChi2PerClusterTPC(4);
2311       SetAcceptKinkDaughters(kFALSE);
2312       SetRequireTPCRefit(kTRUE);
2313       
2314       //ITS
2315       SetRequireITSRefit(kTRUE);
2316       SetMaxDCAToVertexZ(2);
2317       SetMaxDCAToVertexXY();
2318       SetDCAToVertex2D(kFALSE);
2319       
2320       break;
2321     }
2322       
2323     case kLooseCut:
2324     {
2325       AliInfo(Form("Track cuts for matching: Loose cut w/o DCA cut"));
2326       SetMinNClustersTPC(50);
2327       SetAcceptKinkDaughters(kTRUE);
2328       
2329       break;
2330     }
2331   }
2332 }
2333
2334
2335 //________________________________________________________________________
2336 void AliEMCALRecoUtils::SetClusterMatchedToTrack(const AliVEvent *event)
2337 {
2338   // Checks if tracks are matched to EMC clusters and set the matched EMCAL cluster index to ESD track. 
2339
2340   Int_t nTracks = event->GetNumberOfTracks();
2341   for (Int_t iTrack = 0; iTrack < nTracks; ++iTrack) 
2342   {
2343     AliVTrack* track = dynamic_cast<AliVTrack*>(event->GetTrack(iTrack));
2344     if (!track) 
2345     {
2346       AliWarning(Form("Could not receive track %d", iTrack));
2347       continue;
2348     }
2349     
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);
2354     if (esdtrack) { 
2355       if(matchClusIndex != -1) 
2356         esdtrack->SetStatus(AliESDtrack::kEMCALmatch);
2357       else
2358         esdtrack->ResetStatus(AliESDtrack::kEMCALmatch);
2359     } else {
2360       AliAODTrack* aodtrack = dynamic_cast<AliAODTrack*>(track);
2361       if(matchClusIndex != -1) 
2362         aodtrack->SetStatus(AliESDtrack::kEMCALmatch);
2363       else
2364         aodtrack->ResetStatus(AliESDtrack::kEMCALmatch);
2365     }
2366
2367   }
2368   AliDebug(2,"Track matched to closest cluster");  
2369 }
2370
2371 //_________________________________________________________________________
2372 void AliEMCALRecoUtils::SetTracksMatchedToCluster(const AliVEvent *event)
2373 {
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.
2376   
2377   for (Int_t iClus=0; iClus < event->GetNumberOfCaloClusters(); ++iClus) 
2378   {
2379     AliVCluster *cluster = event->GetCaloCluster(iClus);
2380     if (!cluster->IsEMCAL()) 
2381       continue;
2382     
2383     Int_t nTracks = event->GetNumberOfTracks();
2384     TArrayI arrayTrackMatched(nTracks);
2385     
2386     // Get the closest track matched to the cluster
2387     Int_t nMatched = 0;
2388     Int_t matchTrackIndex = GetMatchedTrackIndex(iClus);
2389     if (matchTrackIndex != -1) 
2390     {
2391       arrayTrackMatched[nMatched] = matchTrackIndex;
2392       nMatched++;
2393     }
2394     
2395     // Get all other tracks matched to the cluster
2396     for(Int_t iTrk=0; iTrk<nTracks; ++iTrk) 
2397     {
2398       AliVTrack* track = dynamic_cast<AliVTrack*>(event->GetTrack(iTrk));
2399       if(iTrk == matchTrackIndex) continue;
2400       if(track->GetEMCALcluster() == iClus)
2401       {
2402         arrayTrackMatched[nMatched] = iTrk;
2403         ++nMatched;
2404       }
2405     }
2406     
2407     //printf("Tender::SetTracksMatchedToCluster - cluster E %f, N matches %d, first match %d\n",cluster->E(),nMatched,arrayTrackMatched[0]);
2408     
2409     arrayTrackMatched.Set(nMatched);
2410     AliESDCaloCluster *esdcluster = dynamic_cast<AliESDCaloCluster*>(cluster);
2411     if (esdcluster) 
2412       esdcluster->AddTracksMatched(arrayTrackMatched);
2413     else if (nMatched>0) {
2414       AliAODCaloCluster *aodcluster = dynamic_cast<AliAODCaloCluster*>(cluster);
2415       if (aodcluster)
2416         aodcluster->AddTrackMatched(event->GetTrack(arrayTrackMatched.At(0)));
2417     }
2418     
2419     Float_t eta= -999, phi = -999;
2420     if (matchTrackIndex != -1) 
2421       GetMatchedResiduals(iClus, eta, phi);
2422     cluster->SetTrackDistance(phi, eta);
2423   }
2424   
2425   AliDebug(2,"Cluster matched to tracks");  
2426 }
2427
2428 //___________________________________________________
2429 void AliEMCALRecoUtils::Print(const Option_t *) const 
2430 {
2431   // Print Parameters
2432   
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]);
2440   
2441   printf("Position Recalculation option %d, Particle Type %d, fW0 %2.2f, Recalibrate Data %d \n",fPosAlgo,fParticleType,fW0, fRecalibration);
2442
2443   printf("Matching criteria: ");
2444   if(fCutEtaPhiSum)
2445     {
2446       printf("sqrt(dEta^2+dPhi^2)<%4.3f\n",fCutR);
2447     }
2448   else if(fCutEtaPhiSeparate)
2449     {
2450       printf("dEta<%4.3f, dPhi<%4.3f\n",fCutEta,fCutPhi);
2451     }
2452   else
2453     {
2454       printf("Error\n");
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");
2458     }
2459
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);
2462
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);
2471 }
2472