]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALRecoUtils.cxx
Update from Salvatore to have time recalib working in the tender. Here is the logic:
[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 kPi0MCv3:
680     {
681       //Same as beam test corrected, change parameters
682       //fNonLinearityParams[0] =  9.81039e-01
683       //fNonLinearityParams[1] =  1.13508e-01;
684       //fNonLinearityParams[2] =  1.00173e+00; 
685       //fNonLinearityParams[3] =  9.67998e-02;
686       //fNonLinearityParams[4] =  2.19381e+02;
687       //fNonLinearityParams[5] =  6.31604e+01;
688       //fNonLinearityParams[6] =  1;
689       energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
690       
691       break;
692     }
693       
694       
695     case kPi0GammaGamma:
696     {
697       //Non-Linearity correction (from Olga Data with function p0+p1*exp(-p2*E))
698       //fNonLinearityParams[0] = 1.04;
699       //fNonLinearityParams[1] = -0.1445;
700       //fNonLinearityParams[2] = 1.046;
701       energy /= (fNonLinearityParams[0]+fNonLinearityParams[1]*exp(-fNonLinearityParams[2]*energy)); //Olga function
702       break;
703     }
704       
705     case kPi0GammaConversion:
706     {
707       //Non-Linearity correction (Nicolas from Dimitri Data with function C*[1-a*exp(-b*E)])
708       //fNonLinearityParams[0] = 0.139393/0.1349766;
709       //fNonLinearityParams[1] = 0.0566186;
710       //fNonLinearityParams[2] = 0.982133;
711       energy /= fNonLinearityParams[0]*(1-fNonLinearityParams[1]*exp(-fNonLinearityParams[2]*energy));
712       
713       break;
714     }
715       
716     case kBeamTest:
717     {
718       //From beam test, Alexei's results, for different ZS thresholds
719       //                        th=30 MeV; th = 45 MeV; th = 75 MeV
720       //fNonLinearityParams[0] = 1.007;      1.003;      1.002 
721       //fNonLinearityParams[1] = 0.894;      0.719;      0.797 
722       //fNonLinearityParams[2] = 0.246;      0.334;      0.358 
723       //Rescale the param[0] with 1.03
724       energy /= fNonLinearityParams[0]/(1+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]));
725       
726       break;
727     }
728       
729     case kBeamTestCorrected:
730     {
731       //From beam test, corrected for material between beam and EMCAL
732       //fNonLinearityParams[0] =  0.99078
733       //fNonLinearityParams[1] =  0.161499;
734       //fNonLinearityParams[2] =  0.655166; 
735       //fNonLinearityParams[3] =  0.134101;
736       //fNonLinearityParams[4] =  163.282;
737       //fNonLinearityParams[5] =  23.6904;
738       //fNonLinearityParams[6] =  0.978;
739         energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
740
741       break;
742     }
743       
744     case kNoCorrection:
745       AliDebug(2,"No correction on the energy\n");
746       break;
747       
748   }
749
750   return energy;
751 }
752
753 //__________________________________________________
754 void AliEMCALRecoUtils::InitNonLinearityParam()
755 {
756   //Initialising Non Linearity Parameters
757   
758   if(fNonLinearityFunction == kPi0MC) 
759   {
760     fNonLinearityParams[0] = 1.014;
761     fNonLinearityParams[1] = -0.03329;
762     fNonLinearityParams[2] = -0.3853;
763     fNonLinearityParams[3] = 0.5423;
764     fNonLinearityParams[4] = -0.4335;
765   }
766   
767   if(fNonLinearityFunction == kPi0MCv2) 
768   {
769     fNonLinearityParams[0] = 3.11111e-02;
770     fNonLinearityParams[1] =-5.71666e-02; 
771     fNonLinearityParams[2] = 5.67995e-01;      
772   }
773   
774   if(fNonLinearityFunction == kPi0MCv3) 
775   {
776     fNonLinearityParams[0] =  9.81039e-01;
777     fNonLinearityParams[1] =  1.13508e-01;
778     fNonLinearityParams[2] =  1.00173e+00; 
779     fNonLinearityParams[3] =  9.67998e-02;
780     fNonLinearityParams[4] =  2.19381e+02;
781     fNonLinearityParams[5] =  6.31604e+01;
782     fNonLinearityParams[6] =  1;
783   }
784   
785   if(fNonLinearityFunction == kPi0GammaGamma) 
786   {
787     fNonLinearityParams[0] = 1.04;
788     fNonLinearityParams[1] = -0.1445;
789     fNonLinearityParams[2] = 1.046;
790   }  
791
792   if(fNonLinearityFunction == kPi0GammaConversion) 
793   {
794     fNonLinearityParams[0] = 0.139393;
795     fNonLinearityParams[1] = 0.0566186;
796     fNonLinearityParams[2] = 0.982133;
797   }  
798
799   if(fNonLinearityFunction == kBeamTest) 
800   {
801     if(fNonLinearThreshold == 30) 
802     {
803       fNonLinearityParams[0] = 1.007; 
804       fNonLinearityParams[1] = 0.894; 
805       fNonLinearityParams[2] = 0.246; 
806     }
807     if(fNonLinearThreshold == 45) 
808     {
809       fNonLinearityParams[0] = 1.003; 
810       fNonLinearityParams[1] = 0.719; 
811       fNonLinearityParams[2] = 0.334; 
812     }
813     if(fNonLinearThreshold == 75) 
814     {
815       fNonLinearityParams[0] = 1.002; 
816       fNonLinearityParams[1] = 0.797; 
817       fNonLinearityParams[2] = 0.358; 
818     }
819   }
820
821   if(fNonLinearityFunction == kBeamTestCorrected) 
822   {
823     fNonLinearityParams[0] =  0.99078;
824     fNonLinearityParams[1] =  0.161499;
825     fNonLinearityParams[2] =  0.655166; 
826     fNonLinearityParams[3] =  0.134101;
827     fNonLinearityParams[4] =  163.282;
828     fNonLinearityParams[5] =  23.6904;
829     fNonLinearityParams[6] =  0.978;
830   }
831 }
832
833 //_________________________________________________________
834 Float_t  AliEMCALRecoUtils::GetDepth(const Float_t energy, 
835                                      const Int_t iParticle, 
836                                      const Int_t iSM) const 
837 {
838   //Calculate shower depth for a given cluster energy and particle type
839
840   // parameters 
841   Float_t x0    = 1.31;
842   Float_t ecr   = 8;
843   Float_t depth = 0;
844   
845   switch ( iParticle )
846   {
847     case kPhoton:
848       depth = x0 * (TMath::Log(energy*1000/ ecr) + 0.5); //Multiply energy by 1000 to transform to MeV
849       break;
850       
851     case kElectron:
852       depth = x0 * (TMath::Log(energy*1000/ ecr) - 0.5); //Multiply energy by 1000 to transform to MeV
853       break;
854       
855     case kHadron:
856       // hadron 
857       // boxes anc. here
858       if(gGeoManager)
859       {
860         gGeoManager->cd("ALIC_1/XEN1_1");
861         TGeoNode        *geoXEn1    = gGeoManager->GetCurrentNode();
862         TGeoNodeMatrix  *geoSM      = dynamic_cast<TGeoNodeMatrix *>(geoXEn1->GetDaughter(iSM));
863         if(geoSM)
864         {
865           TGeoVolume      *geoSMVol   = geoSM->GetVolume(); 
866           TGeoShape       *geoSMShape = geoSMVol->GetShape();
867           TGeoBBox        *geoBox     = dynamic_cast<TGeoBBox *>(geoSMShape);
868           if(geoBox) depth = 0.5 * geoBox->GetDX()*2 ;
869           else AliFatal("Null GEANT box");
870         }
871         else AliFatal("NULL  GEANT node matrix");
872       }
873       else
874       {//electron
875         depth = x0 * (TMath::Log(energy*1000 / ecr)  - 0.5); //Multiply energy by 1000 to transform to MeV
876       }
877         
878       break;
879       
880     default://photon
881       depth = x0 * (TMath::Log(energy*1000 / ecr) + 0.5); //Multiply energy by 1000 to transform to MeV
882   }  
883   
884   return depth;
885 }
886
887 //____________________________________________________________________
888 void AliEMCALRecoUtils::GetMaxEnergyCell(const AliEMCALGeometry *geom, 
889                                          AliVCaloCells* cells, 
890                                          const AliVCluster* clu, 
891                                          Int_t  & absId,  
892                                          Int_t  & iSupMod, 
893                                          Int_t  & ieta, 
894                                          Int_t  & iphi, 
895                                          Bool_t & shared)
896 {
897   //For a given CaloCluster gets the absId of the cell 
898   //with maximum energy deposit.
899   
900   Double_t eMax        = -1.;
901   Double_t eCell       = -1.;
902   Float_t  fraction    = 1.;
903   Float_t  recalFactor = 1.;
904   Int_t    cellAbsId   = -1 ;
905
906   Int_t iTower  = -1;
907   Int_t iIphi   = -1;
908   Int_t iIeta   = -1;
909   Int_t iSupMod0= -1;
910
911   if(!clu)
912   {
913     AliInfo("Cluster pointer null!");
914     absId=-1; iSupMod0=-1, ieta = -1; iphi = -1; shared = -1;
915     return;
916   }
917   
918   for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++) 
919   {
920     cellAbsId = clu->GetCellAbsId(iDig);
921     fraction  = clu->GetCellAmplitudeFraction(iDig);
922     //printf("a Cell %d, id, %d, amp %f, fraction %f\n",iDig,cellAbsId,cells->GetCellAmplitude(cellAbsId),fraction);
923     if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
924     geom->GetCellIndex(cellAbsId,iSupMod,iTower,iIphi,iIeta); 
925     geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);
926     if     (iDig==0) 
927     {
928       iSupMod0=iSupMod;
929     }
930     else if(iSupMod0!=iSupMod) 
931     {
932       shared = kTRUE;
933       //printf("AliEMCALRecoUtils::GetMaxEnergyCell() - SHARED CLUSTER\n");
934     }
935     if(!fCellsRecalibrated && IsRecalibrationOn()) 
936     {
937       recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
938     }
939     eCell  = cells->GetCellAmplitude(cellAbsId)*fraction*recalFactor;
940     //printf("b Cell %d, id, %d, amp %f, fraction %f\n",iDig,cellAbsId,eCell,fraction);
941     if(eCell > eMax)  
942     { 
943       eMax  = eCell; 
944       absId = cellAbsId;
945       //printf("\t new max: cell %d, e %f, ecell %f\n",maxId, eMax,eCell);
946     }
947   }// cell loop
948   
949   //Get from the absid the supermodule, tower and eta/phi numbers
950   geom->GetCellIndex(absId,iSupMod,iTower,iIphi,iIeta); 
951   //Gives SuperModule and Tower numbers
952   geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,
953                                          iIphi, iIeta,iphi,ieta); 
954   //printf("Max id %d, iSM %d, col %d, row %d\n",absId,iSupMod,ieta,iphi);
955   //printf("Max end---\n");
956 }
957
958 //______________________________________
959 void AliEMCALRecoUtils::InitParameters()
960 {
961   // Initialize data members with default values
962   
963   fParticleType = kPhoton;
964   fPosAlgo      = kUnchanged;
965   fW0           = 4.5;
966   
967   fNonLinearityFunction = kNoCorrection;
968   fNonLinearThreshold   = 30;
969   
970   fExoticCellFraction     = 0.97;
971   fExoticCellDiffTime     = 1e6;
972   fExoticCellMinAmplitude = 0.5;
973   
974   fAODFilterMask = 32;
975   
976   fCutEtaPhiSum      = kTRUE;
977   fCutEtaPhiSeparate = kFALSE;
978   
979   fCutR   = 0.05; 
980   fCutEta = 0.025; 
981   fCutPhi = 0.05;
982   
983   fClusterWindow = 100;
984   fMass          = 0.139;
985   
986   fStepSurface   = 20.;                      
987   fStepCluster   = 5.;
988   fTrackCutsType = kLooseCut;
989   
990   fCutMinTrackPt     = 0;
991   fCutMinNClusterTPC = -1;
992   fCutMinNClusterITS = -1;
993   
994   fCutMaxChi2PerClusterTPC  = 1e10;
995   fCutMaxChi2PerClusterITS  = 1e10;
996   
997   fCutRequireTPCRefit     = kFALSE;
998   fCutRequireITSRefit     = kFALSE;
999   fCutAcceptKinkDaughters = kFALSE;
1000   
1001   fCutMaxDCAToVertexXY = 1e10;             
1002   fCutMaxDCAToVertexZ  = 1e10;              
1003   fCutDCAToVertex2D    = kFALSE;
1004   
1005   
1006   //Misalignment matrices
1007   for(Int_t i = 0; i < 15 ; i++) 
1008   {
1009     fMisalTransShift[i] = 0.; 
1010     fMisalRotShift[i]   = 0.; 
1011   }
1012   
1013   //Non linearity
1014   for(Int_t i = 0; i < 7  ; i++) fNonLinearityParams[i] = 0.; 
1015   
1016   //For kBeamTestCorrected case, but default is no correction
1017   fNonLinearityParams[0] =  0.99078;
1018   fNonLinearityParams[1] =  0.161499;
1019   fNonLinearityParams[2] =  0.655166; 
1020   fNonLinearityParams[3] =  0.134101;
1021   fNonLinearityParams[4] =  163.282;
1022   fNonLinearityParams[5] =  23.6904;
1023   fNonLinearityParams[6] =  0.978;
1024   
1025   //For kPi0GammaGamma case
1026   //fNonLinearityParams[0] = 0.1457/0.1349766/1.038;
1027   //fNonLinearityParams[1] = -0.02024/0.1349766/1.038;
1028   //fNonLinearityParams[2] = 1.046;
1029   
1030   //Cluster energy smearing
1031   fSmearClusterEnergy   = kFALSE;
1032   fSmearClusterParam[0] = 0.07; // * sqrt E term
1033   fSmearClusterParam[1] = 0.00; // * E term
1034   fSmearClusterParam[2] = 0.00; // constant
1035 }
1036
1037 //_____________________________________________________
1038 void AliEMCALRecoUtils::InitEMCALRecalibrationFactors()
1039 {
1040   //Init EMCAL recalibration factors
1041   AliDebug(2,"AliCalorimeterUtils::InitEMCALRecalibrationFactors()");
1042   //In order to avoid rewriting the same histograms
1043   Bool_t oldStatus = TH1::AddDirectoryStatus();
1044   TH1::AddDirectory(kFALSE);
1045   
1046   fEMCALRecalibrationFactors = new TObjArray(12);
1047   for (int i = 0; i < 12; i++) 
1048     fEMCALRecalibrationFactors->Add(new TH2F(Form("EMCALRecalFactors_SM%d",i),
1049                                              Form("EMCALRecalFactors_SM%d",i),  48, 0, 48, 24, 0, 24));
1050   //Init the histograms with 1
1051   for (Int_t sm = 0; sm < 12; sm++) 
1052   {
1053     for (Int_t i = 0; i < 48; i++) 
1054     {
1055       for (Int_t j = 0; j < 24; j++) 
1056       {
1057         SetEMCALChannelRecalibrationFactor(sm,i,j,1.);
1058       }
1059     }
1060   }
1061   
1062   fEMCALRecalibrationFactors->SetOwner(kTRUE);
1063   fEMCALRecalibrationFactors->Compress();
1064   
1065   //In order to avoid rewriting the same histograms
1066   TH1::AddDirectory(oldStatus);    
1067 }
1068
1069 //_________________________________________________________
1070 void AliEMCALRecoUtils::InitEMCALTimeRecalibrationFactors()
1071 {
1072   //Init EMCAL recalibration factors
1073   AliDebug(2,"AliCalorimeterUtils::InitEMCALRecalibrationFactors()");
1074   //In order to avoid rewriting the same histograms
1075   Bool_t oldStatus = TH1::AddDirectoryStatus();
1076   TH1::AddDirectory(kFALSE);
1077   
1078   fEMCALTimeRecalibrationFactors = new TObjArray(4);
1079   for (int i = 0; i < 4; i++) 
1080     fEMCALTimeRecalibrationFactors->Add(new TH1F(Form("hAllTimeAvBC%d",i),
1081                                                  Form("hAllTimeAvBC%d",i),  
1082                                                  48*24*12,0.,48*24*12)          );
1083   //Init the histograms with 1
1084   for (Int_t bc = 0; bc < 4; bc++) 
1085   {
1086     for (Int_t i = 0; i < 48*24*12; i++) 
1087       SetEMCALChannelTimeRecalibrationFactor(bc,i,0.);
1088   }
1089   
1090   fEMCALTimeRecalibrationFactors->SetOwner(kTRUE);
1091   fEMCALTimeRecalibrationFactors->Compress();
1092   
1093   //In order to avoid rewriting the same histograms
1094   TH1::AddDirectory(oldStatus);    
1095 }
1096
1097 //____________________________________________________
1098 void AliEMCALRecoUtils::InitEMCALBadChannelStatusMap()
1099 {
1100   //Init EMCAL bad channels map
1101   AliDebug(2,"AliEMCALRecoUtils::InitEMCALBadChannelStatusMap()");
1102   //In order to avoid rewriting the same histograms
1103   Bool_t oldStatus = TH1::AddDirectoryStatus();
1104   TH1::AddDirectory(kFALSE);
1105   
1106   fEMCALBadChannelMap = new TObjArray(12);
1107   //TH2F * hTemp = new  TH2I("EMCALBadChannelMap","EMCAL SuperModule bad channel map", 48, 0, 48, 24, 0, 24);
1108   for (int i = 0; i < 12; i++) 
1109   {
1110     fEMCALBadChannelMap->Add(new TH2I(Form("EMCALBadChannelMap_Mod%d",i),Form("EMCALBadChannelMap_Mod%d",i), 48, 0, 48, 24, 0, 24));
1111   }
1112   
1113   fEMCALBadChannelMap->SetOwner(kTRUE);
1114   fEMCALBadChannelMap->Compress();
1115   
1116   //In order to avoid rewriting the same histograms
1117   TH1::AddDirectory(oldStatus);    
1118 }
1119
1120 //____________________________________________________________________________
1121 void AliEMCALRecoUtils::RecalibrateClusterEnergy(const AliEMCALGeometry* geom, 
1122                                                  AliVCluster * cluster, 
1123                                                  AliVCaloCells * cells, 
1124                                                  const Int_t bc)
1125 {
1126   // Recalibrate the cluster energy and Time, considering the recalibration map 
1127   // and the energy of the cells and time that compose the cluster.
1128   // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1129   
1130   if(!cluster)
1131   {
1132     AliInfo("Cluster pointer null!");
1133     return;
1134   }  
1135   
1136   //Get the cluster number of cells and list of absId, check what kind of cluster do we have.
1137   UShort_t * index    = cluster->GetCellsAbsId() ;
1138   Double_t * fraction = cluster->GetCellsAmplitudeFraction() ;
1139   Int_t ncells = cluster->GetNCells();
1140   
1141   //Initialize some used variables
1142   Float_t energy = 0;
1143   Int_t   absId  =-1;
1144   Int_t   icol   =-1, irow =-1, imod=1;
1145   Float_t factor = 1, frac = 0;
1146   Int_t   absIdMax = -1;
1147   Float_t emax     = 0;
1148   
1149   //Loop on the cells, get the cell amplitude and recalibration factor, multiply and and to the new energy
1150   for(Int_t icell = 0; icell < ncells; icell++)
1151   {
1152     absId = index[icell];
1153     frac =  fraction[icell];
1154     if(frac < 1e-5) frac = 1; //in case of EMCAL, this is set as 0 since unfolding is off
1155     
1156     if(!fCellsRecalibrated && IsRecalibrationOn()) 
1157     {
1158       // Energy  
1159       Int_t iTower = -1, iIphi = -1, iIeta = -1; 
1160       geom->GetCellIndex(absId,imod,iTower,iIphi,iIeta); 
1161       if(fEMCALRecalibrationFactors->GetEntries() <= imod) continue;
1162       geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,irow,icol);      
1163       factor = GetEMCALChannelRecalibrationFactor(imod,icol,irow);
1164       
1165       AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - recalibrate cell: module %d, col %d, row %d, cell fraction %f,recalibration factor %f, cell energy %f\n",
1166                       imod,icol,irow,frac,factor,cells->GetCellAmplitude(absId)));
1167       
1168     } 
1169     
1170     energy += cells->GetCellAmplitude(absId)*factor*frac;
1171     
1172     if(emax < cells->GetCellAmplitude(absId)*factor*frac)
1173     {
1174       emax     = cells->GetCellAmplitude(absId)*factor*frac;
1175       absIdMax = absId;
1176     }
1177   }
1178   
1179   AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - Energy before %f, after %f \n",cluster->E(),energy));
1180
1181   cluster->SetE(energy);
1182
1183   // Recalculate time of cluster
1184   Double_t timeorg = cluster->GetTOF();
1185
1186   Double_t time = cells->GetCellTime(absIdMax);
1187   if(!fCellsRecalibrated && IsTimeRecalibrationOn())
1188     RecalibrateCellTime(absIdMax,bc,time);
1189
1190   cluster->SetTOF(time);
1191
1192   AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - Time before %f, after %f \n",timeorg,cluster->GetTOF()));
1193 }
1194
1195 //_____________________________________________________________
1196 void AliEMCALRecoUtils::RecalibrateCells(AliVCaloCells * cells,
1197                                          const Int_t bc)
1198 {
1199   // Recalibrate the cells time and energy, considering the recalibration map and the energy 
1200   // of the cells that compose the cluster.
1201   // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1202
1203   if(!IsRecalibrationOn() && !IsTimeRecalibrationOn()) return;
1204   
1205   if(!cells)
1206   {
1207     AliInfo("Cells pointer null!");
1208     return;
1209   }  
1210   
1211   Short_t  absId  =-1;
1212   Bool_t   accept = kFALSE;
1213   Float_t  ecell  = 0;
1214   Double_t tcell  = 0;
1215   Double_t ecellin = 0;
1216   Double_t tcellin = 0;
1217   Short_t  mclabel = -1;
1218   Double_t efrac = 0;
1219   
1220   Int_t nEMcell  = cells->GetNumberOfCells() ;  
1221   for (Int_t iCell = 0; iCell < nEMcell; iCell++) 
1222   { 
1223     cells->GetCell( iCell, absId, ecellin, tcellin, mclabel, efrac );
1224     
1225     accept = AcceptCalibrateCell(absId, bc, ecell ,tcell ,cells); 
1226     if(!accept) 
1227     {
1228       ecell = 0;
1229       tcell = -1;
1230     }
1231     
1232     //Set new values
1233     cells->SetCell(iCell,absId,ecell, tcell, mclabel, efrac);
1234   }
1235
1236   fCellsRecalibrated = kTRUE;
1237 }
1238
1239 //_______________________________________________________________________________________________________
1240 void AliEMCALRecoUtils::RecalibrateCellTime(const Int_t absId, const Int_t bc, Double_t & celltime) const
1241 {
1242   // Recalibrate time of cell with absID  considering the recalibration map 
1243   // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1244   
1245   if(!fCellsRecalibrated && IsTimeRecalibrationOn() && bc >= 0)
1246   {
1247     celltime -= GetEMCALChannelTimeRecalibrationFactor(bc%4,absId)*1.e-9;    ;  
1248   }
1249 }
1250   
1251 //______________________________________________________________________________
1252 void AliEMCALRecoUtils::RecalculateClusterPosition(const AliEMCALGeometry *geom, 
1253                                                    AliVCaloCells* cells, 
1254                                                    AliVCluster* clu)
1255 {
1256   //For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1257   
1258   if(!clu)
1259   {
1260     AliInfo("Cluster pointer null!");
1261     return;
1262   }
1263     
1264   if     (fPosAlgo==kPosTowerGlobal) RecalculateClusterPositionFromTowerGlobal( geom, cells, clu);
1265   else if(fPosAlgo==kPosTowerIndex)  RecalculateClusterPositionFromTowerIndex ( geom, cells, clu);
1266   else   AliDebug(2,"Algorithm to recalculate position not selected, do nothing.");
1267 }  
1268
1269 //_____________________________________________________________________________________________
1270 void AliEMCALRecoUtils::RecalculateClusterPositionFromTowerGlobal(const AliEMCALGeometry *geom, 
1271                                                                   AliVCaloCells* cells, 
1272                                                                   AliVCluster* clu)
1273 {
1274   // For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1275   // The algorithm is a copy of what is done in AliEMCALRecPoint
1276   
1277   Double_t eCell       = 0.;
1278   Float_t  fraction    = 1.;
1279   Float_t  recalFactor = 1.;
1280   
1281   Int_t    absId   = -1;
1282   Int_t    iTower  = -1, iIphi  = -1, iIeta  = -1;
1283   Int_t    iSupModMax = -1, iSM=-1, iphi   = -1, ieta   = -1;
1284   Float_t  weight = 0.,  totalWeight=0.;
1285   Float_t  newPos[3] = {0,0,0};
1286   Double_t pLocal[3], pGlobal[3];
1287   Bool_t shared = kFALSE;
1288
1289   Float_t  clEnergy = clu->E(); //Energy already recalibrated previously
1290   GetMaxEnergyCell(geom, cells, clu, absId,  iSupModMax, ieta, iphi,shared);
1291   Double_t depth = GetDepth(clEnergy,fParticleType,iSupModMax) ;
1292   
1293   //printf("** Cluster energy %f, ncells %d, depth %f\n",clEnergy,clu->GetNCells(),depth);
1294   
1295   for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++) 
1296   {
1297     absId = clu->GetCellAbsId(iDig);
1298     fraction  = clu->GetCellAmplitudeFraction(iDig);
1299     if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1300     
1301     if (!fCellsRecalibrated)
1302     {
1303       geom->GetCellIndex(absId,iSM,iTower,iIphi,iIeta); 
1304       geom->GetCellPhiEtaIndexInSModule(iSM,iTower,iIphi, iIeta,iphi,ieta);      
1305       
1306       if(IsRecalibrationOn()) 
1307       {
1308         recalFactor = GetEMCALChannelRecalibrationFactor(iSM,ieta,iphi);
1309       }
1310     }
1311     
1312     eCell  = cells->GetCellAmplitude(absId)*fraction*recalFactor;
1313     
1314     weight = GetCellWeight(eCell,clEnergy);
1315     totalWeight += weight;
1316     
1317     geom->RelPosCellInSModule(absId,depth,pLocal[0],pLocal[1],pLocal[2]);
1318     //printf("pLocal (%f,%f,%f), SM %d, absId %d\n",pLocal[0],pLocal[1],pLocal[2],iSupModMax,absId);
1319     geom->GetGlobal(pLocal,pGlobal,iSupModMax);
1320     //printf("pLocal (%f,%f,%f)\n",pGlobal[0],pGlobal[1],pGlobal[2]);
1321
1322     for(int i=0; i<3; i++ ) newPos[i] += (weight*pGlobal[i]);
1323   }// cell loop
1324   
1325   if(totalWeight>0)
1326   {
1327     for(int i=0; i<3; i++ )    newPos[i] /= totalWeight;
1328   }
1329     
1330   //Float_t pos[]={0,0,0};
1331   //clu->GetPosition(pos);
1332   //printf("OldPos  : %2.3f,%2.3f,%2.3f\n",pos[0],pos[1],pos[2]);
1333   //printf("NewPos  : %2.3f,%2.3f,%2.3f\n",newPos[0],newPos[1],newPos[2]);
1334   
1335   if(iSupModMax > 1) //sector 1
1336   {
1337     newPos[0] +=fMisalTransShift[3];//-=3.093; 
1338     newPos[1] +=fMisalTransShift[4];//+=6.82;
1339     newPos[2] +=fMisalTransShift[5];//+=1.635;
1340     //printf("   +    : %2.3f,%2.3f,%2.3f\n",fMisalTransShift[3],fMisalTransShift[4],fMisalTransShift[5]);
1341   } else //sector 0
1342   {
1343     newPos[0] +=fMisalTransShift[0];//+=1.134;
1344     newPos[1] +=fMisalTransShift[1];//+=8.2;
1345     newPos[2] +=fMisalTransShift[2];//+=1.197;
1346     //printf("   +    : %2.3f,%2.3f,%2.3f\n",fMisalTransShift[0],fMisalTransShift[1],fMisalTransShift[2]);
1347   }
1348   //printf("NewPos : %2.3f,%2.3f,%2.3f\n",newPos[0],newPos[1],newPos[2]);
1349
1350   clu->SetPosition(newPos);
1351 }  
1352
1353 //____________________________________________________________________________________________
1354 void AliEMCALRecoUtils::RecalculateClusterPositionFromTowerIndex(const AliEMCALGeometry *geom, 
1355                                                                  AliVCaloCells* cells, 
1356                                                                  AliVCluster* clu)
1357 {
1358   // For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1359   // The algorithm works with the tower indeces, averages the indeces and from them it calculates the global position
1360   
1361   Double_t eCell       = 1.;
1362   Float_t  fraction    = 1.;
1363   Float_t  recalFactor = 1.;
1364   
1365   Int_t absId   = -1;
1366   Int_t iTower  = -1;
1367   Int_t iIphi   = -1, iIeta   = -1;
1368   Int_t iSupMod = -1, iSupModMax = -1;
1369   Int_t iphi = -1, ieta =-1;
1370   Bool_t shared = kFALSE;
1371
1372   Float_t clEnergy = clu->E(); //Energy already recalibrated previously.
1373   GetMaxEnergyCell(geom, cells, clu, absId,  iSupModMax, ieta, iphi,shared);
1374   Float_t  depth = GetDepth(clEnergy,fParticleType,iSupMod) ;
1375
1376   Float_t weight = 0., weightedCol = 0., weightedRow = 0., totalWeight=0.;
1377   Bool_t areInSameSM = kTRUE; //exclude clusters with cells in different SMs for now
1378   Int_t startingSM = -1;
1379   
1380   for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++) 
1381   {
1382     absId = clu->GetCellAbsId(iDig);
1383     fraction  = clu->GetCellAmplitudeFraction(iDig);
1384     if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1385
1386     if     (iDig==0)  startingSM = iSupMod;
1387     else if(iSupMod != startingSM) areInSameSM = kFALSE;
1388
1389     eCell  = cells->GetCellAmplitude(absId);
1390     
1391     geom->GetCellIndex(absId,iSupMod,iTower,iIphi,iIeta); 
1392     geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);    
1393     
1394     if (!fCellsRecalibrated)
1395     {
1396       if(IsRecalibrationOn()) 
1397       {
1398         recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1399       }
1400     }
1401     
1402     eCell  = cells->GetCellAmplitude(absId)*fraction*recalFactor;
1403     
1404     weight = GetCellWeight(eCell,clEnergy);
1405     if(weight < 0) weight = 0;
1406     totalWeight += weight;
1407     weightedCol += ieta*weight;
1408     weightedRow += iphi*weight;
1409     
1410     //printf("Max cell? cell %d, amplitude org %f, fraction %f, recalibration %f, amplitude new %f \n",cellAbsId, cells->GetCellAmplitude(cellAbsId), fraction, recalFactor, eCell) ;
1411   }// cell loop
1412     
1413   Float_t xyzNew[]={0.,0.,0.};
1414   if(areInSameSM == kTRUE) 
1415   {
1416     //printf("In Same SM\n");
1417     weightedCol = weightedCol/totalWeight;
1418     weightedRow = weightedRow/totalWeight;
1419     geom->RecalculateTowerPosition(weightedRow, weightedCol, iSupModMax, depth, fMisalTransShift, fMisalRotShift, xyzNew); 
1420   } 
1421   else 
1422   {
1423     //printf("In Different SM\n");
1424     geom->RecalculateTowerPosition(iphi,        ieta,        iSupModMax, depth, fMisalTransShift, fMisalRotShift, xyzNew); 
1425   }
1426   
1427   clu->SetPosition(xyzNew);
1428 }
1429
1430 //___________________________________________________________________________________________
1431 void AliEMCALRecoUtils::RecalculateClusterDistanceToBadChannel(const AliEMCALGeometry * geom, 
1432                                                                AliVCaloCells* cells, 
1433                                                                AliVCluster * cluster)
1434 {           
1435   //re-evaluate distance to bad channel with updated bad map
1436   
1437   if(!fRecalDistToBadChannels) return;
1438   
1439   if(!cluster)
1440   {
1441     AliInfo("Cluster pointer null!");
1442     return;
1443   }  
1444   
1445   //Get channels map of the supermodule where the cluster is.
1446   Int_t absIdMax  = -1, iSupMod =-1, icolM = -1, irowM = -1;
1447   Bool_t shared = kFALSE;
1448   GetMaxEnergyCell(geom, cells, cluster, absIdMax,  iSupMod, icolM, irowM, shared);
1449   TH2D* hMap  = (TH2D*)fEMCALBadChannelMap->At(iSupMod);
1450
1451   Int_t dRrow, dRcol;  
1452   Float_t  minDist = 10000.;
1453   Float_t  dist    = 0.;
1454   
1455   //Loop on tower status map 
1456   for(Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++)
1457   {
1458     for(Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++)
1459     {
1460       //Check if tower is bad.
1461       if(hMap->GetBinContent(icol,irow)==0) continue;
1462       //printf("AliEMCALRecoUtils::RecalculateDistanceToBadChannels() - \n \t Bad channel in SM %d, col %d, row %d, \n \t Cluster max in col %d, row %d\n",
1463       //       iSupMod,icol, irow, icolM,irowM);
1464       
1465       dRrow=TMath::Abs(irowM-irow);
1466       dRcol=TMath::Abs(icolM-icol);
1467       dist=TMath::Sqrt(dRrow*dRrow+dRcol*dRcol);
1468       if(dist < minDist)
1469       {
1470         //printf("MIN DISTANCE TO BAD %2.2f\n",dist);
1471         minDist = dist;
1472       }
1473     }
1474   }
1475   
1476   //In case the cluster is shared by 2 SuperModules, need to check the map of the second Super Module
1477   if (shared) 
1478   {
1479     TH2D* hMap2 = 0;
1480     Int_t iSupMod2 = -1;
1481     
1482     //The only possible combinations are (0,1), (2,3) ... (8,9)
1483     if(iSupMod%2) iSupMod2 = iSupMod-1;
1484     else          iSupMod2 = iSupMod+1;
1485     hMap2  = (TH2D*)fEMCALBadChannelMap->At(iSupMod2);
1486     
1487     //Loop on tower status map of second super module
1488     for(Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++)
1489     {
1490       for(Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++)
1491       {
1492         //Check if tower is bad.
1493         if(hMap2->GetBinContent(icol,irow)==0) continue;
1494         //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",
1495         //     iSupMod2,icol, irow,iSupMod,icolM,irowM);
1496         dRrow=TMath::Abs(irow-irowM);
1497         
1498         if(iSupMod%2) 
1499         {
1500           dRcol=TMath::Abs(icol-(AliEMCALGeoParams::fgkEMCALCols+icolM));
1501         } else 
1502         {
1503           dRcol=TMath::Abs(AliEMCALGeoParams::fgkEMCALCols+icol-icolM);
1504         }                    
1505         
1506         dist=TMath::Sqrt(dRrow*dRrow+dRcol*dRcol);
1507         if(dist < minDist) minDist = dist;        
1508       }
1509     }
1510   }// shared cluster in 2 SuperModules
1511   
1512   AliDebug(2,Form("Max cluster cell (SM,col,row)=(%d %d %d) - Distance to Bad Channel %2.2f",iSupMod, icolM, irowM, minDist));
1513   cluster->SetDistanceToBadChannel(minDist);
1514 }
1515
1516 //__________________________________________________________________
1517 void AliEMCALRecoUtils::RecalculateClusterPID(AliVCluster * cluster)
1518 {           
1519   //re-evaluate identification parameters with bayesian
1520   
1521   if(!cluster)
1522   {
1523     AliInfo("Cluster pointer null!");
1524     return;
1525   }
1526   
1527   if ( cluster->GetM02() != 0)
1528     fPIDUtils->ComputePID(cluster->E(),cluster->GetM02());
1529   
1530   Float_t pidlist[AliPID::kSPECIESCN+1];
1531   for(Int_t i = 0; i < AliPID::kSPECIESCN+1; i++) pidlist[i] = fPIDUtils->GetPIDFinal(i);
1532         
1533   cluster->SetPID(pidlist);
1534 }
1535
1536 //___________________________________________________________________________________________________________________
1537 void AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(const AliEMCALGeometry * geom, 
1538                                                                 AliVCaloCells* cells, 
1539                                                                 AliVCluster * cluster,
1540                                                                 Float_t & l0,   Float_t & l1,   
1541                                                                 Float_t & disp, Float_t & dEta, Float_t & dPhi,
1542                                                                 Float_t & sEta, Float_t & sPhi, Float_t & sEtaPhi)
1543 {
1544   // Calculates new center of gravity in the local EMCAL-module coordinates 
1545   // and tranfers into global ALICE coordinates
1546   // Calculates Dispersion and main axis
1547   
1548   if(!cluster)
1549   {
1550     AliInfo("Cluster pointer null!");
1551     return;
1552   }
1553     
1554   Double_t eCell       = 0.;
1555   Float_t  fraction    = 1.;
1556   Float_t  recalFactor = 1.;
1557
1558   Int_t    iSupMod = -1;
1559   Int_t    iTower  = -1;
1560   Int_t    iIphi   = -1;
1561   Int_t    iIeta   = -1;
1562   Int_t    iphi    = -1;
1563   Int_t    ieta    = -1;
1564   Double_t etai    = -1.;
1565   Double_t phii    = -1.;
1566   
1567   Int_t    nstat   = 0 ;
1568   Float_t  wtot    = 0.;
1569   Double_t w       = 0.;
1570   Double_t etaMean = 0.;
1571   Double_t phiMean = 0.;
1572     
1573   //Loop on cells
1574   for(Int_t iDigit=0; iDigit < cluster->GetNCells(); iDigit++) 
1575   {
1576     //Get from the absid the supermodule, tower and eta/phi numbers
1577     geom->GetCellIndex(cluster->GetCellAbsId(iDigit),iSupMod,iTower,iIphi,iIeta); 
1578     geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi,iIeta, iphi,ieta);        
1579     
1580     //Get the cell energy, if recalibration is on, apply factors
1581     fraction  = cluster->GetCellAmplitudeFraction(iDigit);
1582     if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1583     
1584     if (!fCellsRecalibrated)
1585     {
1586       if(IsRecalibrationOn()) 
1587       {
1588         recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1589       }
1590     }
1591     
1592     eCell  = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1593     
1594     if(cluster->E() > 0 && eCell > 0)
1595     {
1596       w  = GetCellWeight(eCell,cluster->E());
1597       
1598       etai=(Double_t)ieta;
1599       phii=(Double_t)iphi;  
1600       
1601       if(w > 0.0) 
1602       {
1603         wtot += w ;
1604         nstat++;            
1605         //Shower shape
1606         sEta     += w * etai * etai ;
1607         etaMean  += w * etai ;
1608         sPhi     += w * phii * phii ;
1609         phiMean  += w * phii ; 
1610         sEtaPhi  += w * etai * phii ; 
1611       }
1612     } 
1613     else
1614       AliError(Form("Wrong energy %f and/or amplitude %f\n", eCell, cluster->E()));
1615   }//cell loop
1616   
1617   //Normalize to the weight  
1618   if (wtot > 0) 
1619   {
1620     etaMean /= wtot ;
1621     phiMean /= wtot ;
1622   }
1623   else
1624     AliError(Form("Wrong weight %f\n", wtot));
1625   
1626   //Calculate dispersion  
1627   for(Int_t iDigit=0; iDigit < cluster->GetNCells(); iDigit++) 
1628   {
1629     //Get from the absid the supermodule, tower and eta/phi numbers
1630     geom->GetCellIndex(cluster->GetCellAbsId(iDigit),iSupMod,iTower,iIphi,iIeta); 
1631     geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi,iIeta, iphi,ieta);
1632     
1633     //Get the cell energy, if recalibration is on, apply factors
1634     fraction  = cluster->GetCellAmplitudeFraction(iDigit);
1635     if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1636     if (IsRecalibrationOn()) 
1637     {
1638       recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1639     }
1640     eCell  = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1641     
1642     if(cluster->E() > 0 && eCell > 0)
1643     {
1644       w  = GetCellWeight(eCell,cluster->E());
1645       
1646       etai=(Double_t)ieta;
1647       phii=(Double_t)iphi;    
1648       if(w > 0.0) 
1649       { 
1650         disp +=  w *((etai-etaMean)*(etai-etaMean)+(phii-phiMean)*(phii-phiMean)); 
1651         dEta +=  w * (etai-etaMean)*(etai-etaMean) ; 
1652         dPhi +=  w * (phii-phiMean)*(phii-phiMean) ; 
1653       }
1654     }
1655     else
1656       AliError(Form("Wrong energy %f and/or amplitude %f\n", eCell, cluster->E()));
1657   }// cell loop
1658   
1659   //Normalize to the weigth and set shower shape parameters
1660   if (wtot > 0 && nstat > 1) 
1661   {
1662     disp    /= wtot ;
1663     dEta    /= wtot ;
1664     dPhi    /= wtot ;
1665     sEta    /= wtot ;
1666     sPhi    /= wtot ;
1667     sEtaPhi /= wtot ;
1668     
1669     sEta    -= etaMean * etaMean ;
1670     sPhi    -= phiMean * phiMean ;
1671     sEtaPhi -= etaMean * phiMean ;
1672     
1673     l0 = (0.5 * (sEta + sPhi) + TMath::Sqrt( 0.25 * (sEta - sPhi) * (sEta - sPhi) + sEtaPhi * sEtaPhi ));
1674     l1 = (0.5 * (sEta + sPhi) - TMath::Sqrt( 0.25 * (sEta - sPhi) * (sEta - sPhi) + sEtaPhi * sEtaPhi ));
1675   }
1676   else
1677   {
1678     l0   = 0. ;
1679     l1   = 0. ;
1680     dEta = 0. ; dPhi = 0. ; disp    = 0. ;
1681     sEta = 0. ; sPhi = 0. ; sEtaPhi = 0. ;
1682   }  
1683   
1684 }
1685
1686 //____________________________________________________________________________________________
1687 void AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(const AliEMCALGeometry * geom, 
1688                                                                 AliVCaloCells* cells, 
1689                                                                 AliVCluster * cluster)
1690 {
1691   // Calculates new center of gravity in the local EMCAL-module coordinates 
1692   // and tranfers into global ALICE coordinates
1693   // Calculates Dispersion and main axis and puts them into the cluster
1694   
1695   Float_t l0   = 0., l1   = 0.;
1696   Float_t disp = 0., dEta = 0., dPhi    = 0.; 
1697   Float_t sEta = 0., sPhi = 0., sEtaPhi = 0.;
1698   
1699   AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(geom,cells,cluster,l0,l1,disp,
1700                                                              dEta, dPhi, sEta, sPhi, sEtaPhi);
1701   
1702   cluster->SetM02(l0);
1703   cluster->SetM20(l1);
1704   if(disp > 0. ) cluster->SetDispersion(TMath::Sqrt(disp)) ;
1705   
1706
1707
1708 //____________________________________________________________________________
1709 void AliEMCALRecoUtils::FindMatches(AliVEvent *event,
1710                                     TObjArray * clusterArr,  
1711                                     const AliEMCALGeometry *geom)
1712 {
1713   //This function should be called before the cluster loop
1714   //Before call this function, please recalculate the cluster positions
1715   //Given the input event, loop over all the tracks, select the closest cluster as matched with fCutR
1716   //Store matched cluster indexes and residuals
1717   
1718   fMatchedTrackIndex  ->Reset();
1719   fMatchedClusterIndex->Reset();
1720   fResidualPhi->Reset();
1721   fResidualEta->Reset();
1722   
1723   fMatchedTrackIndex  ->Set(1000);
1724   fMatchedClusterIndex->Set(1000);
1725   fResidualPhi->Set(1000);
1726   fResidualEta->Set(1000);
1727   
1728   AliESDEvent* esdevent = dynamic_cast<AliESDEvent*> (event);
1729   AliAODEvent* aodevent = dynamic_cast<AliAODEvent*> (event);
1730   
1731   // init the magnetic field if not already on
1732   if(!TGeoGlobalMagField::Instance()->GetField())
1733   {
1734     AliInfo("Init the magnetic field\n");
1735     if     (esdevent) 
1736     {
1737       esdevent->InitMagneticField();
1738     }
1739     else if(aodevent)
1740     {
1741       Double_t curSol = 30000*aodevent->GetMagneticField()/5.00668;
1742       Double_t curDip = 6000 *aodevent->GetMuonMagFieldScale();
1743       AliMagF *field  = AliMagF::CreateFieldMap(curSol,curDip);
1744       TGeoGlobalMagField::Instance()->SetField(field);
1745     }
1746     else
1747     {
1748       AliInfo("Mag Field not initialized, null esd/aod evetn pointers");
1749     }
1750     
1751   } // Init mag field
1752   
1753   TObjArray *clusterArray = 0x0;
1754   if(!clusterArr)
1755     {
1756       clusterArray = new TObjArray(event->GetNumberOfCaloClusters());
1757       for(Int_t icl=0; icl<event->GetNumberOfCaloClusters(); icl++)
1758   {
1759     AliVCluster *cluster = (AliVCluster*) event->GetCaloCluster(icl);
1760     if(geom && !IsGoodCluster(cluster,geom,(AliVCaloCells*)event->GetEMCALCells())) continue;
1761     clusterArray->AddAt(cluster,icl);
1762   }
1763     }
1764   
1765   Int_t    matched=0;
1766   Double_t cv[21];
1767   for (Int_t i=0; i<21;i++) cv[i]=0;
1768   for(Int_t itr=0; itr<event->GetNumberOfTracks(); itr++)
1769   {
1770     AliExternalTrackParam *trackParam = 0;
1771     
1772     //If the input event is ESD, the starting point for extrapolation is TPCOut, if available, or TPCInner 
1773     AliESDtrack *esdTrack = 0;
1774     AliAODTrack *aodTrack = 0;
1775     if(esdevent)
1776     {
1777       esdTrack = esdevent->GetTrack(itr);
1778       if(!esdTrack) continue;
1779       if(!IsAccepted(esdTrack)) continue;
1780       if(esdTrack->Pt()<fCutMinTrackPt) continue;
1781       Double_t phi = esdTrack->Phi()*TMath::RadToDeg();
1782       if(TMath::Abs(esdTrack->Eta())>0.8 || phi <= 20 || phi >= 240 ) continue;
1783       trackParam =  const_cast<AliExternalTrackParam*>(esdTrack->GetInnerParam());
1784     }
1785     
1786     //If the input event is AOD, the starting point for extrapolation is at vertex
1787     //AOD tracks are selected according to its filterbit.
1788     else if(aodevent)
1789     {
1790       aodTrack = aodevent->GetTrack(itr);
1791       if(!aodTrack) continue;
1792       if(!aodTrack->TestFilterMask(fAODFilterMask)) continue; //Select AOD tracks that fulfill GetStandardITSTPCTrackCuts2010()
1793       if(aodTrack->Pt()<fCutMinTrackPt) continue;
1794       Double_t phi = aodTrack->Phi()*TMath::RadToDeg();
1795       if(TMath::Abs(aodTrack->Eta())>0.8 || phi <= 20 || phi >= 240 ) continue;
1796       Double_t pos[3],mom[3];
1797       aodTrack->GetXYZ(pos);
1798       aodTrack->GetPxPyPz(mom);
1799       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()));
1800       trackParam= new AliExternalTrackParam(pos,mom,cv,aodTrack->Charge());
1801     }
1802     
1803     //Return if the input data is not "AOD" or "ESD"
1804     else
1805     {
1806       printf("Wrong input data type! Should be \"AOD\" or \"ESD\"\n");
1807       if(clusterArray)
1808   {
1809     clusterArray->Clear();
1810     delete clusterArray;
1811   }
1812       return;
1813     }
1814     
1815     if(!trackParam) continue;
1816
1817     //Extrapolate the track to EMCal surface
1818     AliExternalTrackParam emcalParam(*trackParam);
1819     Float_t eta, phi;
1820     if(!ExtrapolateTrackToEMCalSurface(&emcalParam, 430., fMass, fStepSurface, eta, phi)) 
1821       {
1822   if(aodevent && trackParam) delete trackParam;
1823   continue;
1824       }
1825
1826 //    if(esdevent)
1827 //      {
1828 //  esdTrack->SetOuterParam(&emcalParam,AliExternalTrackParam::kMultSec);
1829 //      }
1830
1831     if(TMath::Abs(eta)>0.75 || (phi) < 70*TMath::DegToRad() || (phi) > 190*TMath::DegToRad())
1832       {
1833   if(aodevent && trackParam) delete trackParam;
1834   continue;
1835       }
1836
1837
1838     //Find matched clusters
1839     Int_t index = -1;
1840     Float_t dEta = -999, dPhi = -999;
1841     if(!clusterArr)
1842       {
1843   index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArray, dEta, dPhi);  
1844       }
1845     else
1846       {
1847   index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArr, dEta, dPhi);  
1848       }  
1849     
1850     if(index>-1)
1851     {
1852       fMatchedTrackIndex   ->AddAt(itr,matched);
1853       fMatchedClusterIndex ->AddAt(index,matched);
1854       fResidualEta         ->AddAt(dEta,matched);
1855       fResidualPhi         ->AddAt(dPhi,matched);
1856       matched++;
1857     }
1858     if(aodevent && trackParam) delete trackParam;
1859   }//track loop
1860
1861   if(clusterArray)
1862     {
1863       clusterArray->Clear();
1864       delete clusterArray;
1865     }
1866   
1867   AliDebug(2,Form("Number of matched pairs = %d !\n",matched));
1868   
1869   fMatchedTrackIndex   ->Set(matched);
1870   fMatchedClusterIndex ->Set(matched);
1871   fResidualPhi         ->Set(matched);
1872   fResidualEta         ->Set(matched);
1873 }
1874
1875 //________________________________________________________________________________
1876 Int_t AliEMCALRecoUtils::FindMatchedClusterInEvent(const AliESDtrack *track, 
1877                                                    const AliVEvent *event, 
1878                                                    const AliEMCALGeometry *geom, 
1879                                                    Float_t &dEta, Float_t &dPhi)
1880 {
1881   //
1882   // This function returns the index of matched cluster to input track
1883   // Returns -1 if no match is found
1884   Int_t index = -1;
1885   Double_t phiV = track->Phi()*TMath::RadToDeg();
1886   if(TMath::Abs(track->Eta())>0.8 || phiV <= 20 || phiV >= 240 ) return index;
1887   AliExternalTrackParam *trackParam = const_cast<AliExternalTrackParam*>(track->GetInnerParam());
1888   if(!trackParam) return index;
1889   AliExternalTrackParam emcalParam(*trackParam);
1890   Float_t eta, phi;
1891   if(!ExtrapolateTrackToEMCalSurface(&emcalParam, 430., fMass, fStepSurface, eta, phi)) return index;
1892   if(TMath::Abs(eta)>0.75 || (phi) < 70*TMath::DegToRad() || (phi) > 190*TMath::DegToRad()) return index;
1893
1894   TObjArray *clusterArr = new TObjArray(event->GetNumberOfCaloClusters());
1895
1896   for(Int_t icl=0; icl<event->GetNumberOfCaloClusters(); icl++)
1897   {
1898     AliVCluster *cluster = (AliVCluster*) event->GetCaloCluster(icl);
1899     if(geom && !IsGoodCluster(cluster,geom,(AliVCaloCells*)event->GetEMCALCells())) continue;
1900     clusterArr->AddAt(cluster,icl);
1901   }
1902
1903   index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArr, dEta, dPhi);  
1904   clusterArr->Clear();
1905   delete clusterArr;
1906   
1907   return index;
1908 }
1909
1910 //_______________________________________________________________________________________________
1911 Int_t  AliEMCALRecoUtils::FindMatchedClusterInClusterArr(const AliExternalTrackParam *emcalParam, 
1912                                                          AliExternalTrackParam *trkParam, 
1913                                                          const TObjArray * clusterArr, 
1914                                                          Float_t &dEta, Float_t &dPhi)
1915 {
1916   // Find matched cluster in array
1917   
1918   dEta=-999, dPhi=-999;
1919   Float_t dRMax = fCutR, dEtaMax=fCutEta, dPhiMax=fCutPhi;
1920   Int_t index = -1;
1921   Float_t tmpEta=-999, tmpPhi=-999;
1922
1923   Double_t exPos[3] = {0.,0.,0.};
1924   if(!emcalParam->GetXYZ(exPos)) return index;
1925
1926   Float_t clsPos[3] = {0.,0.,0.};
1927   for(Int_t icl=0; icl<clusterArr->GetEntriesFast(); icl++)
1928     {
1929       AliVCluster *cluster = dynamic_cast<AliVCluster*> (clusterArr->At(icl)) ;
1930       if(!cluster || !cluster->IsEMCAL()) continue;
1931       cluster->GetPosition(clsPos);
1932       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));
1933       if(dR > fClusterWindow) continue;
1934
1935       AliExternalTrackParam trkPamTmp (*trkParam);//Retrieve the starting point every time before the extrapolation
1936       if(!ExtrapolateTrackToCluster(&trkPamTmp, cluster, fMass, fStepCluster, tmpEta, tmpPhi)) continue;
1937       if(fCutEtaPhiSum)
1938         {
1939           Float_t tmpR=TMath::Sqrt(tmpEta*tmpEta + tmpPhi*tmpPhi);
1940           if(tmpR<dRMax)
1941       {
1942         dRMax=tmpR;
1943         dEtaMax=tmpEta;
1944         dPhiMax=tmpPhi;
1945         index=icl;
1946       }
1947         }
1948       else if(fCutEtaPhiSeparate)
1949         {
1950           if(TMath::Abs(tmpEta)<TMath::Abs(dEtaMax) && TMath::Abs(tmpPhi)<TMath::Abs(dPhiMax))
1951       {
1952         dEtaMax = tmpEta;
1953         dPhiMax = tmpPhi;
1954         index=icl;
1955       }
1956         }
1957       else
1958         {
1959           printf("Error: please specify your cut criteria\n");
1960           printf("To cut on sqrt(dEta^2+dPhi^2), use: SwitchOnCutEtaPhiSum()\n");
1961           printf("To cut on dEta and dPhi separately, use: SwitchOnCutEtaPhiSeparate()\n");
1962           return index;
1963         }
1964     }
1965
1966   dEta=dEtaMax;
1967   dPhi=dPhiMax;
1968
1969   return index;
1970 }
1971
1972 //------------------------------------------------------------------------------------
1973 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToEMCalSurface(AliExternalTrackParam *trkParam, 
1974                                                          const Double_t emcalR,
1975                                                          const Double_t mass, 
1976                                                          const Double_t step, 
1977                                                          Float_t &eta, 
1978                                                          Float_t &phi)
1979 {
1980   //Extrapolate track to EMCAL surface
1981   
1982   eta = -999, phi = -999;
1983   if(!trkParam) return kFALSE;
1984   if(!AliTrackerBase::PropagateTrackToBxByBz(trkParam, emcalR, mass, step, kTRUE, 0.8, -1)) return kFALSE;
1985   Double_t trkPos[3] = {0.,0.,0.};
1986   if(!trkParam->GetXYZ(trkPos)) return kFALSE;
1987   TVector3 trkPosVec(trkPos[0],trkPos[1],trkPos[2]);
1988   eta = trkPosVec.Eta();
1989   phi = trkPosVec.Phi();
1990   if(phi<0)
1991     phi += 2*TMath::Pi();
1992
1993   return kTRUE;
1994 }
1995
1996 //-----------------------------------------------------------------------------------
1997 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToPosition(AliExternalTrackParam *trkParam, 
1998                                                      const Float_t *clsPos, 
1999                                                      Double_t mass, 
2000                                                      Double_t step, 
2001                                                      Float_t &tmpEta, 
2002                                                      Float_t &tmpPhi)
2003 {
2004   //
2005   //Return the residual by extrapolating a track param to a global position
2006   //
2007   tmpEta = -999;
2008   tmpPhi = -999;
2009   if(!trkParam) return kFALSE;
2010   Double_t trkPos[3] = {0.,0.,0.};
2011   TVector3 vec(clsPos[0],clsPos[1],clsPos[2]);
2012   Double_t alpha =  ((int)(vec.Phi()*TMath::RadToDeg()/20)+0.5)*20*TMath::DegToRad();
2013   vec.RotateZ(-alpha); //Rotate the cluster to the local extrapolation coordinate system
2014   if(!AliTrackerBase::PropagateTrackToBxByBz(trkParam, vec.X(), mass, step,kTRUE, 0.8, -1)) return kFALSE;
2015   if(!trkParam->GetXYZ(trkPos)) return kFALSE; //Get the extrapolated global position
2016
2017   TVector3 clsPosVec(clsPos[0],clsPos[1],clsPos[2]);
2018   TVector3 trkPosVec(trkPos[0],trkPos[1],trkPos[2]);
2019
2020   // track cluster matching
2021   tmpPhi = clsPosVec.DeltaPhi(trkPosVec);    // tmpPhi is between -pi and pi
2022   tmpEta = clsPosVec.Eta()-trkPosVec.Eta();
2023
2024   return kTRUE;
2025 }
2026
2027 //----------------------------------------------------------------------------------
2028 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToCluster(AliExternalTrackParam *trkParam, 
2029                                                     const AliVCluster *cluster, 
2030                                                     const Double_t mass, 
2031                                                     const Double_t step, 
2032                                                     Float_t &tmpEta, 
2033                                                     Float_t &tmpPhi)
2034 {
2035   //
2036   //Return the residual by extrapolating a track param to a cluster
2037   //
2038   tmpEta = -999;
2039   tmpPhi = -999;
2040   if(!cluster || !trkParam) return kFALSE;
2041
2042   Float_t clsPos[3] = {0.,0.,0.};
2043   cluster->GetPosition(clsPos);
2044
2045   return ExtrapolateTrackToPosition(trkParam, clsPos, mass, step, tmpEta, tmpPhi);
2046 }
2047
2048 //---------------------------------------------------------------------------------
2049 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToCluster(AliExternalTrackParam *trkParam, 
2050                                                     const AliVCluster *cluster, 
2051                                                     Float_t &tmpEta, 
2052                                                     Float_t &tmpPhi)
2053 {
2054   //
2055   //Return the residual by extrapolating a track param to a clusterfStepCluster
2056   //
2057
2058   return ExtrapolateTrackToCluster(trkParam, cluster, fMass, fStepCluster, tmpEta, tmpPhi);
2059 }
2060
2061 //_______________________________________________________________________
2062 void AliEMCALRecoUtils::GetMatchedResiduals(const Int_t clsIndex, 
2063                                             Float_t &dEta, Float_t &dPhi)
2064 {
2065   //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2066   //Get the residuals dEta and dPhi for this cluster to the closest track
2067   //Works with ESDs and AODs
2068
2069   if( FindMatchedPosForCluster(clsIndex) >= 999 )
2070   {
2071     AliDebug(2,"No matched tracks found!\n");
2072     dEta=999.;
2073     dPhi=999.;
2074     return;
2075   }
2076   dEta = fResidualEta->At(FindMatchedPosForCluster(clsIndex));
2077   dPhi = fResidualPhi->At(FindMatchedPosForCluster(clsIndex));
2078 }
2079
2080 //______________________________________________________________________________________________
2081 void AliEMCALRecoUtils::GetMatchedClusterResiduals(Int_t trkIndex, Float_t &dEta, Float_t &dPhi)
2082 {
2083   //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2084   //Get the residuals dEta and dPhi for this track to the closest cluster
2085   //Works with ESDs and AODs
2086
2087   if( FindMatchedPosForTrack(trkIndex) >= 999 )
2088   {
2089     AliDebug(2,"No matched cluster found!\n");
2090     dEta=999.;
2091     dPhi=999.;
2092     return;
2093   }
2094   dEta = fResidualEta->At(FindMatchedPosForTrack(trkIndex));
2095   dPhi = fResidualPhi->At(FindMatchedPosForTrack(trkIndex));
2096 }
2097
2098 //__________________________________________________________
2099 Int_t AliEMCALRecoUtils::GetMatchedTrackIndex(Int_t clsIndex)
2100 {
2101   //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2102   //Get the index of matched track to this cluster
2103   //Works with ESDs and AODs
2104   
2105   if(IsClusterMatched(clsIndex))
2106     return fMatchedTrackIndex->At(FindMatchedPosForCluster(clsIndex));
2107   else 
2108     return -1; 
2109 }
2110
2111 //__________________________________________________________
2112 Int_t AliEMCALRecoUtils::GetMatchedClusterIndex(Int_t trkIndex)
2113 {
2114   //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2115   //Get the index of matched cluster to this track
2116   //Works with ESDs and AODs
2117   
2118   if(IsTrackMatched(trkIndex))
2119     return fMatchedClusterIndex->At(FindMatchedPosForTrack(trkIndex));
2120   else 
2121     return -1; 
2122 }
2123
2124 //______________________________________________________________
2125 Bool_t AliEMCALRecoUtils::IsClusterMatched(Int_t clsIndex) const
2126 {
2127   //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2128   //Returns if the cluster has a match
2129   if(FindMatchedPosForCluster(clsIndex) < 999) 
2130     return kTRUE;
2131   else
2132     return kFALSE;
2133 }
2134
2135 //____________________________________________________________
2136 Bool_t AliEMCALRecoUtils::IsTrackMatched(Int_t trkIndex) const 
2137 {
2138   //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2139   //Returns if the track has a match
2140   if(FindMatchedPosForTrack(trkIndex) < 999) 
2141     return kTRUE;
2142   else
2143     return kFALSE;
2144 }
2145
2146 //______________________________________________________________________
2147 UInt_t AliEMCALRecoUtils::FindMatchedPosForCluster(Int_t clsIndex) const
2148 {
2149   //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2150   //Returns the position of the match in the fMatchedClusterIndex array
2151   Float_t tmpR = fCutR;
2152   UInt_t pos = 999;
2153   
2154   for(Int_t i=0; i<fMatchedClusterIndex->GetSize(); i++) 
2155   {
2156     if(fMatchedClusterIndex->At(i)==clsIndex) 
2157     {
2158       Float_t r = TMath::Sqrt(fResidualEta->At(i)*fResidualEta->At(i) + fResidualPhi->At(i)*fResidualPhi->At(i));
2159       if(r<tmpR) 
2160       {
2161         pos=i;
2162         tmpR=r;
2163         AliDebug(3,Form("Matched cluster index: index: %d, dEta: %2.4f, dPhi: %2.4f.\n",
2164                         fMatchedClusterIndex->At(i),fResidualEta->At(i),fResidualPhi->At(i)));
2165       }
2166     }
2167   }
2168   return pos;
2169 }
2170
2171 //____________________________________________________________________
2172 UInt_t AliEMCALRecoUtils::FindMatchedPosForTrack(Int_t trkIndex) const
2173 {
2174   //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2175   //Returns the position of the match in the fMatchedTrackIndex array
2176   Float_t tmpR = fCutR;
2177   UInt_t pos = 999;
2178   
2179   for(Int_t i=0; i<fMatchedTrackIndex->GetSize(); i++) 
2180   {
2181     if(fMatchedTrackIndex->At(i)==trkIndex) 
2182     {
2183       Float_t r = TMath::Sqrt(fResidualEta->At(i)*fResidualEta->At(i) + fResidualPhi->At(i)*fResidualPhi->At(i));
2184       if(r<tmpR) 
2185       {
2186         pos=i;
2187         tmpR=r;
2188         AliDebug(3,Form("Matched track index: index: %d, dEta: %2.4f, dPhi: %2.4f.\n",
2189                         fMatchedTrackIndex->At(i),fResidualEta->At(i),fResidualPhi->At(i)));
2190       }
2191     }
2192   }
2193   return pos;
2194 }
2195
2196 //__________________________________________________________________________
2197 Bool_t AliEMCALRecoUtils::IsGoodCluster(AliVCluster *cluster, 
2198                                         const AliEMCALGeometry *geom, 
2199                                         AliVCaloCells* cells,const Int_t bc)
2200 {
2201   // check if the cluster survives some quality cut
2202   //
2203   //
2204   Bool_t isGood=kTRUE;
2205
2206   if(!cluster || !cluster->IsEMCAL())              return kFALSE;
2207   
2208   if(ClusterContainsBadChannel(geom,cluster->GetCellsAbsId(),cluster->GetNCells())) return kFALSE;
2209   
2210   if(!CheckCellFiducialRegion(geom,cluster,cells)) return kFALSE;
2211   
2212   if(IsExoticCluster(cluster, cells,bc))           return kFALSE;
2213
2214   return isGood;
2215 }
2216
2217 //__________________________________________________________
2218 Bool_t AliEMCALRecoUtils::IsAccepted(AliESDtrack *esdTrack)
2219 {
2220   // Given a esd track, return whether the track survive all the cuts
2221
2222   // The different quality parameter are first
2223   // retrieved from the track. then it is found out what cuts the
2224   // track did not survive and finally the cuts are imposed.
2225
2226   UInt_t status = esdTrack->GetStatus();
2227
2228   Int_t nClustersITS = esdTrack->GetITSclusters(0);
2229   Int_t nClustersTPC = esdTrack->GetTPCclusters(0);
2230
2231   Float_t chi2PerClusterITS = -1;
2232   Float_t chi2PerClusterTPC = -1;
2233   if (nClustersITS!=0)
2234     chi2PerClusterITS = esdTrack->GetITSchi2()/Float_t(nClustersITS);
2235   if (nClustersTPC!=0) 
2236     chi2PerClusterTPC = esdTrack->GetTPCchi2()/Float_t(nClustersTPC);
2237
2238
2239   //DCA cuts
2240   if(fTrackCutsType==kGlobalCut)
2241     {
2242       Float_t maxDCAToVertexXYPtDep = 0.0182 + 0.0350/TMath::Power(esdTrack->Pt(),1.01); //This expression comes from AliESDtrackCuts::GetStandardITSTPCTrackCuts2010()
2243       //AliDebug(3,Form("Track pT = %f, DCAtoVertexXY = %f",esdTrack->Pt(),MaxDCAToVertexXYPtDep));
2244       SetMaxDCAToVertexXY(maxDCAToVertexXYPtDep); //Set pT dependent DCA cut to vertex in x-y plane
2245     }
2246
2247
2248   Float_t b[2];
2249   Float_t bCov[3];
2250   esdTrack->GetImpactParameters(b,bCov);
2251   if (bCov[0]<=0 || bCov[2]<=0) 
2252   {
2253     AliDebug(1, "Estimated b resolution lower or equal zero!");
2254     bCov[0]=0; bCov[2]=0;
2255   }
2256
2257   Float_t dcaToVertexXY = b[0];
2258   Float_t dcaToVertexZ = b[1];
2259   Float_t dcaToVertex = -1;
2260
2261   if (fCutDCAToVertex2D)
2262     dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY/fCutMaxDCAToVertexXY/fCutMaxDCAToVertexXY + dcaToVertexZ*dcaToVertexZ/fCutMaxDCAToVertexZ/fCutMaxDCAToVertexZ);
2263   else
2264     dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY + dcaToVertexZ*dcaToVertexZ);
2265     
2266   // cut the track?
2267   
2268   Bool_t cuts[kNCuts];
2269   for (Int_t i=0; i<kNCuts; i++) cuts[i]=kFALSE;
2270   
2271   // track quality cuts
2272   if (fCutRequireTPCRefit && (status&AliESDtrack::kTPCrefit)==0)
2273     cuts[0]=kTRUE;
2274   if (fCutRequireITSRefit && (status&AliESDtrack::kITSrefit)==0)
2275     cuts[1]=kTRUE;
2276   if (nClustersTPC<fCutMinNClusterTPC)
2277     cuts[2]=kTRUE;
2278   if (nClustersITS<fCutMinNClusterITS) 
2279     cuts[3]=kTRUE;
2280   if (chi2PerClusterTPC>fCutMaxChi2PerClusterTPC) 
2281     cuts[4]=kTRUE; 
2282   if (chi2PerClusterITS>fCutMaxChi2PerClusterITS) 
2283     cuts[5]=kTRUE;  
2284   if (!fCutAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0)
2285     cuts[6]=kTRUE;
2286   if (fCutDCAToVertex2D && dcaToVertex > 1)
2287     cuts[7] = kTRUE;
2288   if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexXY) > fCutMaxDCAToVertexXY)
2289     cuts[8] = kTRUE;
2290   if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexZ) > fCutMaxDCAToVertexZ)
2291     cuts[9] = kTRUE;
2292
2293   if(fTrackCutsType==kGlobalCut)
2294     {
2295       //Require at least one SPD point + anything else in ITS
2296       if( (esdTrack->HasPointOnITSLayer(0) || esdTrack->HasPointOnITSLayer(1)) == kFALSE)
2297   cuts[10] = kTRUE;
2298     }
2299
2300   Bool_t cut=kFALSE;
2301   for (Int_t i=0; i<kNCuts; i++)
2302     if (cuts[i]) { cut = kTRUE ; }
2303
2304     // cut the track
2305   if (cut) 
2306     return kFALSE;
2307   else 
2308     return kTRUE;
2309 }
2310
2311 //_____________________________________
2312 void AliEMCALRecoUtils::InitTrackCuts()
2313 {
2314   //Intilize the track cut criteria
2315   //By default these cuts are set according to AliESDtrackCuts::GetStandardTPCOnlyTrackCuts()
2316   //Also you can customize the cuts using the setters
2317   
2318   switch (fTrackCutsType)
2319   {
2320     case kTPCOnlyCut:
2321     {
2322       AliInfo(Form("Track cuts for matching: GetStandardTPCOnlyTrackCuts()"));
2323       //TPC
2324       SetMinNClustersTPC(70);
2325       SetMaxChi2PerClusterTPC(4);
2326       SetAcceptKinkDaughters(kFALSE);
2327       SetRequireTPCRefit(kFALSE);
2328       
2329       //ITS
2330       SetRequireITSRefit(kFALSE);
2331       SetMaxDCAToVertexZ(3.2);
2332       SetMaxDCAToVertexXY(2.4);
2333       SetDCAToVertex2D(kTRUE);
2334       
2335       break;
2336     }
2337       
2338     case kGlobalCut:
2339     {
2340       AliInfo(Form("Track cuts for matching: GetStandardITSTPCTrackCuts2010(kTURE)"));
2341       //TPC
2342       SetMinNClustersTPC(70);
2343       SetMaxChi2PerClusterTPC(4);
2344       SetAcceptKinkDaughters(kFALSE);
2345       SetRequireTPCRefit(kTRUE);
2346       
2347       //ITS
2348       SetRequireITSRefit(kTRUE);
2349       SetMaxDCAToVertexZ(2);
2350       SetMaxDCAToVertexXY();
2351       SetDCAToVertex2D(kFALSE);
2352       
2353       break;
2354     }
2355       
2356     case kLooseCut:
2357     {
2358       AliInfo(Form("Track cuts for matching: Loose cut w/o DCA cut"));
2359       SetMinNClustersTPC(50);
2360       SetAcceptKinkDaughters(kTRUE);
2361       
2362       break;
2363     }
2364   }
2365 }
2366
2367
2368 //________________________________________________________________________
2369 void AliEMCALRecoUtils::SetClusterMatchedToTrack(const AliVEvent *event)
2370 {
2371   // Checks if tracks are matched to EMC clusters and set the matched EMCAL cluster index to ESD track. 
2372
2373   Int_t nTracks = event->GetNumberOfTracks();
2374   for (Int_t iTrack = 0; iTrack < nTracks; ++iTrack) 
2375   {
2376     AliVTrack* track = dynamic_cast<AliVTrack*>(event->GetTrack(iTrack));
2377     if (!track) 
2378     {
2379       AliWarning(Form("Could not receive track %d", iTrack));
2380       continue;
2381     }
2382     
2383     Int_t matchClusIndex = GetMatchedClusterIndex(iTrack);       
2384     track->SetEMCALcluster(matchClusIndex); //sets -1 if track not matched within residual
2385     /*the following can be done better if AliVTrack::SetStatus will be there. Patch pending with Andreas/Peter*/
2386     AliESDtrack* esdtrack = dynamic_cast<AliESDtrack*>(track);
2387     if (esdtrack) { 
2388       if(matchClusIndex != -1) 
2389         esdtrack->SetStatus(AliESDtrack::kEMCALmatch);
2390       else
2391         esdtrack->ResetStatus(AliESDtrack::kEMCALmatch);
2392     } else {
2393       AliAODTrack* aodtrack = dynamic_cast<AliAODTrack*>(track);
2394       if(matchClusIndex != -1) 
2395         aodtrack->SetStatus(AliESDtrack::kEMCALmatch);
2396       else
2397         aodtrack->ResetStatus(AliESDtrack::kEMCALmatch);
2398     }
2399
2400   }
2401   AliDebug(2,"Track matched to closest cluster");  
2402 }
2403
2404 //_________________________________________________________________________
2405 void AliEMCALRecoUtils::SetTracksMatchedToCluster(const AliVEvent *event)
2406 {
2407   // Checks if EMC clusters are matched to ESD track.
2408   // Adds track indexes of all the tracks matched to a cluster withing residuals in ESDCalocluster.
2409   
2410   for (Int_t iClus=0; iClus < event->GetNumberOfCaloClusters(); ++iClus) 
2411   {
2412     AliVCluster *cluster = event->GetCaloCluster(iClus);
2413     if (!cluster->IsEMCAL()) 
2414       continue;
2415     
2416     Int_t nTracks = event->GetNumberOfTracks();
2417     TArrayI arrayTrackMatched(nTracks);
2418     
2419     // Get the closest track matched to the cluster
2420     Int_t nMatched = 0;
2421     Int_t matchTrackIndex = GetMatchedTrackIndex(iClus);
2422     if (matchTrackIndex != -1) 
2423     {
2424       arrayTrackMatched[nMatched] = matchTrackIndex;
2425       nMatched++;
2426     }
2427     
2428     // Get all other tracks matched to the cluster
2429     for(Int_t iTrk=0; iTrk<nTracks; ++iTrk) 
2430     {
2431       AliVTrack* track = dynamic_cast<AliVTrack*>(event->GetTrack(iTrk));
2432       if(iTrk == matchTrackIndex) continue;
2433       if(track->GetEMCALcluster() == iClus)
2434       {
2435         arrayTrackMatched[nMatched] = iTrk;
2436         ++nMatched;
2437       }
2438     }
2439     
2440     //printf("Tender::SetTracksMatchedToCluster - cluster E %f, N matches %d, first match %d\n",cluster->E(),nMatched,arrayTrackMatched[0]);
2441     
2442     arrayTrackMatched.Set(nMatched);
2443     AliESDCaloCluster *esdcluster = dynamic_cast<AliESDCaloCluster*>(cluster);
2444     if (esdcluster) 
2445       esdcluster->AddTracksMatched(arrayTrackMatched);
2446     else if (nMatched>0) {
2447       AliAODCaloCluster *aodcluster = dynamic_cast<AliAODCaloCluster*>(cluster);
2448       if (aodcluster)
2449         aodcluster->AddTrackMatched(event->GetTrack(arrayTrackMatched.At(0)));
2450     }
2451     
2452     Float_t eta= -999, phi = -999;
2453     if (matchTrackIndex != -1) 
2454       GetMatchedResiduals(iClus, eta, phi);
2455     cluster->SetTrackDistance(phi, eta);
2456   }
2457   
2458   AliDebug(2,"Cluster matched to tracks");  
2459 }
2460
2461 //___________________________________________________
2462 void AliEMCALRecoUtils::Print(const Option_t *) const 
2463 {
2464   // Print Parameters
2465   
2466   printf("AliEMCALRecoUtils Settings: \n");
2467   printf("Misalignment shifts\n");
2468   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, 
2469                                   fMisalTransShift[i*3],fMisalTransShift[i*3+1],fMisalTransShift[i*3+2],
2470                                   fMisalRotShift[i*3],  fMisalRotShift[i*3+1],  fMisalRotShift[i*3+2]   );
2471   printf("Non linearity function %d, parameters:\n", fNonLinearityFunction);
2472   for(Int_t i=0; i<6; i++) printf("param[%d]=%f\n",i, fNonLinearityParams[i]);
2473   
2474   printf("Position Recalculation option %d, Particle Type %d, fW0 %2.2f, Recalibrate Data %d \n",fPosAlgo,fParticleType,fW0, fRecalibration);
2475
2476   printf("Matching criteria: ");
2477   if(fCutEtaPhiSum)
2478     {
2479       printf("sqrt(dEta^2+dPhi^2)<%4.3f\n",fCutR);
2480     }
2481   else if(fCutEtaPhiSeparate)
2482     {
2483       printf("dEta<%4.3f, dPhi<%4.3f\n",fCutEta,fCutPhi);
2484     }
2485   else
2486     {
2487       printf("Error\n");
2488       printf("please specify your cut criteria\n");
2489       printf("To cut on sqrt(dEta^2+dPhi^2), use: SwitchOnCutEtaPhiSum()\n");
2490       printf("To cut on dEta and dPhi separately, use: SwitchOnCutEtaPhiSeparate()\n");
2491     }
2492
2493   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);
2494   printf("Cluster selection window: dR < %2.0f\n",fClusterWindow);
2495
2496   printf("Track cuts: \n");
2497   printf("Minimum track pT: %1.2f\n",fCutMinTrackPt);
2498   printf("AOD track selection mask: %d\n",fAODFilterMask);
2499   printf("TPCRefit = %d, ITSRefit = %d\n",fCutRequireTPCRefit,fCutRequireITSRefit);
2500   printf("AcceptKinks = %d\n",fCutAcceptKinkDaughters);
2501   printf("MinNCulsterTPC = %d, MinNClusterITS = %d\n",fCutMinNClusterTPC,fCutMinNClusterITS);
2502   printf("MaxChi2TPC = %2.2f, MaxChi2ITS = %2.2f\n",fCutMaxChi2PerClusterTPC,fCutMaxChi2PerClusterITS);
2503   printf("DCSToVertex2D = %d, MaxDCAToVertexXY = %2.2f, MaxDCAToVertexZ = %2.2f\n",fCutDCAToVertex2D,fCutMaxDCAToVertexXY,fCutMaxDCAToVertexZ);
2504 }
2505