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