]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALRecoUtils.cxx
Preparation to do trackmatching on AOD with tender
[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 "AliEMCALCalibTimeDepCorrection.h" // Run dependent
60 #include "AliEMCALPIDUtils.h"
61
62
63 ClassImp(AliEMCALRecoUtils)
64   
65 //_____________________________________
66 AliEMCALRecoUtils::AliEMCALRecoUtils():
67   fParticleType(0),                       fPosAlgo(0),                            fW0(0), 
68   fNonLinearityFunction(0),               fNonLinearThreshold(0),
69   fSmearClusterEnergy(kFALSE),            fRandom(),
70   fCellsRecalibrated(kFALSE),             fRecalibration(kFALSE),                 fEMCALRecalibrationFactors(),
71   fTimeRecalibration(kFALSE),             fEMCALTimeRecalibrationFactors(),
72   fUseRunCorrectionFactors(kFALSE),       fRunCorrectionFactorsSet(kFALSE),
73   fRemoveBadChannels(kFALSE),             fRecalDistToBadChannels(kFALSE),        fEMCALBadChannelMap(),
74   fNCellsFromEMCALBorder(0),              fNoEMCALBorderAtEta0(kTRUE),
75   fRejectExoticCluster(kFALSE),           fRejectExoticCells(kFALSE), 
76   fExoticCellFraction(0),                 fExoticCellDiffTime(0),                 fExoticCellMinAmplitude(0),
77   fPIDUtils(),                            fAODFilterMask(0),
78   fMatchedTrackIndex(0x0),                fMatchedClusterIndex(0x0), 
79   fResidualEta(0x0), fResidualPhi(0x0),   fCutEtaPhiSum(kFALSE),                  fCutEtaPhiSeparate(kFALSE), 
80   fCutR(0),                               fCutEta(0),                             fCutPhi(0),
81   fClusterWindow(0),                      fMass(0),                           
82   fStepSurface(0),                        fStepCluster(0),
83   fTrackCutsType(0),                      fCutMinTrackPt(0),                      fCutMinNClusterTPC(0), 
84   fCutMinNClusterITS(0),                  fCutMaxChi2PerClusterTPC(0),            fCutMaxChi2PerClusterITS(0),
85   fCutRequireTPCRefit(kFALSE),            fCutRequireITSRefit(kFALSE),            fCutAcceptKinkDaughters(kFALSE),
86   fCutMaxDCAToVertexXY(0),                fCutMaxDCAToVertexZ(0),                 fCutDCAToVertex2D(kFALSE)
87 {
88 //
89   // Constructor.
90   // Initialize all constant values which have to be used
91   // during Reco algorithm execution
92   //
93   
94   // Init parameters
95   InitParameters();
96   
97   //Track matching
98   fMatchedTrackIndex     = new TArrayI();
99   fMatchedClusterIndex   = new TArrayI();
100   fResidualPhi           = new TArrayF();
101   fResidualEta           = new TArrayF();
102   fPIDUtils              = new AliEMCALPIDUtils();
103
104   InitTrackCuts();
105 }
106
107 //______________________________________________________________________
108 AliEMCALRecoUtils::AliEMCALRecoUtils(const AliEMCALRecoUtils & reco) 
109 : TNamed(reco), 
110   fParticleType(reco.fParticleType),                         fPosAlgo(reco.fPosAlgo),     fW0(reco.fW0),
111   fNonLinearityFunction(reco.fNonLinearityFunction),         fNonLinearThreshold(reco.fNonLinearThreshold),
112   fSmearClusterEnergy(reco.fSmearClusterEnergy),             fRandom(),
113   fCellsRecalibrated(reco.fCellsRecalibrated),
114   fRecalibration(reco.fRecalibration),                       fEMCALRecalibrationFactors(reco.fEMCALRecalibrationFactors),
115   fTimeRecalibration(reco.fTimeRecalibration),               fEMCALTimeRecalibrationFactors(reco.fEMCALTimeRecalibrationFactors),
116   fUseRunCorrectionFactors(reco.fUseRunCorrectionFactors),   fRunCorrectionFactorsSet(reco.fRunCorrectionFactorsSet),
117   fRemoveBadChannels(reco.fRemoveBadChannels),               fRecalDistToBadChannels(reco.fRecalDistToBadChannels),
118   fEMCALBadChannelMap(reco.fEMCALBadChannelMap),
119   fNCellsFromEMCALBorder(reco.fNCellsFromEMCALBorder),       fNoEMCALBorderAtEta0(reco.fNoEMCALBorderAtEta0),
120   fRejectExoticCluster(reco.fRejectExoticCluster),           fRejectExoticCells(reco.fRejectExoticCells), 
121   fExoticCellFraction(reco.fExoticCellFraction),             fExoticCellDiffTime(reco.fExoticCellDiffTime),               
122   fExoticCellMinAmplitude(reco.fExoticCellMinAmplitude),
123   fPIDUtils(reco.fPIDUtils),                                 fAODFilterMask(reco.fAODFilterMask),
124   fMatchedTrackIndex(  reco.fMatchedTrackIndex?  new TArrayI(*reco.fMatchedTrackIndex):0x0),
125   fMatchedClusterIndex(reco.fMatchedClusterIndex?new TArrayI(*reco.fMatchedClusterIndex):0x0),
126   fResidualEta(        reco.fResidualEta?        new TArrayF(*reco.fResidualEta):0x0),
127   fResidualPhi(        reco.fResidualPhi?        new TArrayF(*reco.fResidualPhi):0x0),
128   fCutEtaPhiSum(reco.fCutEtaPhiSum),                         fCutEtaPhiSeparate(reco.fCutEtaPhiSeparate), 
129   fCutR(reco.fCutR),        fCutEta(reco.fCutEta),           fCutPhi(reco.fCutPhi),
130   fClusterWindow(reco.fClusterWindow),
131   fMass(reco.fMass),        fStepSurface(reco.fStepSurface), fStepCluster(reco.fStepCluster),
132   fTrackCutsType(reco.fTrackCutsType),                       fCutMinTrackPt(reco.fCutMinTrackPt), 
133   fCutMinNClusterTPC(reco.fCutMinNClusterTPC),               fCutMinNClusterITS(reco.fCutMinNClusterITS), 
134   fCutMaxChi2PerClusterTPC(reco.fCutMaxChi2PerClusterTPC),   fCutMaxChi2PerClusterITS(reco.fCutMaxChi2PerClusterITS),
135   fCutRequireTPCRefit(reco.fCutRequireTPCRefit),             fCutRequireITSRefit(reco.fCutRequireITSRefit),
136   fCutAcceptKinkDaughters(reco.fCutAcceptKinkDaughters),     fCutMaxDCAToVertexXY(reco.fCutMaxDCAToVertexXY),    
137   fCutMaxDCAToVertexZ(reco.fCutMaxDCAToVertexZ),             fCutDCAToVertex2D(reco.fCutDCAToVertex2D)
138 {
139   //Copy ctor
140   
141   for(Int_t i = 0; i < 15 ; i++) { fMisalRotShift[i]      = reco.fMisalRotShift[i]      ; 
142                                    fMisalTransShift[i]    = reco.fMisalTransShift[i]    ; } 
143   for(Int_t i = 0; i < 7  ; i++) { fNonLinearityParams[i] = reco.fNonLinearityParams[i] ; }
144   for(Int_t i = 0; i < 3  ; i++) { fSmearClusterParam[i]  = reco.fSmearClusterParam[i]  ; }
145
146 }
147
148
149 //______________________________________________________________________
150 AliEMCALRecoUtils & AliEMCALRecoUtils::operator = (const AliEMCALRecoUtils & reco) 
151 {
152   //Assignment operator
153   
154   if(this == &reco)return *this;
155   ((TNamed *)this)->operator=(reco);
156
157   for(Int_t i = 0; i < 15 ; i++) { fMisalTransShift[i]    = reco.fMisalTransShift[i]    ; 
158                                    fMisalRotShift[i]      = reco.fMisalRotShift[i]      ; }
159   for(Int_t i = 0; i < 7  ; i++) { fNonLinearityParams[i] = reco.fNonLinearityParams[i] ; }
160   for(Int_t i = 0; i < 3  ; i++) { fSmearClusterParam[i]  = reco.fSmearClusterParam[i]  ; }   
161   
162   fParticleType              = reco.fParticleType;
163   fPosAlgo                   = reco.fPosAlgo; 
164   fW0                        = reco.fW0;
165   
166   fNonLinearityFunction      = reco.fNonLinearityFunction;
167   fNonLinearThreshold        = reco.fNonLinearThreshold;
168   fSmearClusterEnergy        = reco.fSmearClusterEnergy;
169
170   fCellsRecalibrated         = reco.fCellsRecalibrated;
171   fRecalibration             = reco.fRecalibration;
172   fEMCALRecalibrationFactors = reco.fEMCALRecalibrationFactors;
173
174   fTimeRecalibration             = reco.fTimeRecalibration;
175   fEMCALTimeRecalibrationFactors = reco.fEMCALTimeRecalibrationFactors;
176
177   fUseRunCorrectionFactors   = reco.fUseRunCorrectionFactors;
178   fRunCorrectionFactorsSet   = reco.fRunCorrectionFactorsSet;
179   
180   fRemoveBadChannels         = reco.fRemoveBadChannels;
181   fRecalDistToBadChannels    = reco.fRecalDistToBadChannels;
182   fEMCALBadChannelMap        = reco.fEMCALBadChannelMap;
183   
184   fNCellsFromEMCALBorder     = reco.fNCellsFromEMCALBorder;
185   fNoEMCALBorderAtEta0       = reco.fNoEMCALBorderAtEta0;
186   
187   fRejectExoticCluster       = reco.fRejectExoticCluster;           
188   fRejectExoticCells         = reco.fRejectExoticCells; 
189   fExoticCellFraction        = reco.fExoticCellFraction;
190   fExoticCellDiffTime        = reco.fExoticCellDiffTime;              
191   fExoticCellMinAmplitude    = reco.fExoticCellMinAmplitude;
192   
193   fPIDUtils                  = reco.fPIDUtils;
194
195   fAODFilterMask             = reco.fAODFilterMask;
196   
197   fCutEtaPhiSum              = reco.fCutEtaPhiSum;
198   fCutEtaPhiSeparate         = reco.fCutEtaPhiSeparate;
199   fCutR                      = reco.fCutR;
200   fCutEta                    = reco.fCutEta;
201   fCutPhi                    = reco.fCutPhi;
202   fClusterWindow             = reco.fClusterWindow;
203   fMass                      = reco.fMass;
204   fStepSurface               = reco.fStepSurface;
205   fStepCluster               = reco.fStepCluster;
206
207   fTrackCutsType             = reco.fTrackCutsType;
208   fCutMinTrackPt             = reco.fCutMinTrackPt;
209   fCutMinNClusterTPC         = reco.fCutMinNClusterTPC;
210   fCutMinNClusterITS         = reco.fCutMinNClusterITS; 
211   fCutMaxChi2PerClusterTPC   = reco.fCutMaxChi2PerClusterTPC;
212   fCutMaxChi2PerClusterITS   = reco.fCutMaxChi2PerClusterITS;
213   fCutRequireTPCRefit        = reco.fCutRequireTPCRefit;
214   fCutRequireITSRefit        = reco.fCutRequireITSRefit;
215   fCutAcceptKinkDaughters    = reco.fCutAcceptKinkDaughters;
216   fCutMaxDCAToVertexXY       = reco.fCutMaxDCAToVertexXY;
217   fCutMaxDCAToVertexZ        = reco.fCutMaxDCAToVertexZ;
218   fCutDCAToVertex2D          = reco.fCutDCAToVertex2D;
219   
220   if(reco.fResidualEta)
221   {
222     // assign or copy construct
223     if(fResidualEta)
224     { 
225       *fResidualEta = *reco.fResidualEta;
226     }
227     else 
228     {
229       fResidualEta = new TArrayF(*reco.fResidualEta);
230     }
231   }
232   else
233   {
234     if(fResidualEta)delete fResidualEta;
235     fResidualEta = 0;
236   }
237   
238   if(reco.fResidualPhi)
239   {
240     // assign or copy construct
241     if(fResidualPhi)
242     { 
243       *fResidualPhi = *reco.fResidualPhi;
244     }
245     else 
246     {
247       fResidualPhi = new TArrayF(*reco.fResidualPhi);
248     }
249   }
250   else
251   {
252     if(fResidualPhi)delete fResidualPhi;
253     fResidualPhi = 0;
254   }
255   
256   if(reco.fMatchedTrackIndex)
257   {
258     // assign or copy construct
259     if(fMatchedTrackIndex)
260     { 
261       *fMatchedTrackIndex = *reco.fMatchedTrackIndex;
262     }
263     else 
264     { 
265       fMatchedTrackIndex = new TArrayI(*reco.fMatchedTrackIndex);
266     }
267   }
268   else
269   {
270     if(fMatchedTrackIndex)delete fMatchedTrackIndex;
271     fMatchedTrackIndex = 0;
272   }  
273   
274   if(reco.fMatchedClusterIndex)
275   {
276     // assign or copy construct
277     if(fMatchedClusterIndex)
278     { 
279       *fMatchedClusterIndex = *reco.fMatchedClusterIndex;
280     }
281     else 
282     {
283       fMatchedClusterIndex = new TArrayI(*reco.fMatchedClusterIndex);
284     }
285   }
286   else
287   {
288     if(fMatchedClusterIndex)delete fMatchedClusterIndex;
289     fMatchedClusterIndex = 0;
290   }
291    
292   return *this;
293 }
294
295
296 //_____________________________________
297 AliEMCALRecoUtils::~AliEMCALRecoUtils()
298 {
299   //Destructor.
300         
301   if(fEMCALRecalibrationFactors) 
302   { 
303     fEMCALRecalibrationFactors->Clear();
304     delete fEMCALRecalibrationFactors;
305   }     
306   
307   if(fEMCALTimeRecalibrationFactors) 
308   { 
309     fEMCALTimeRecalibrationFactors->Clear();
310     delete fEMCALTimeRecalibrationFactors;
311   }     
312   
313   if(fEMCALBadChannelMap) 
314   { 
315     fEMCALBadChannelMap->Clear();
316     delete fEMCALBadChannelMap;
317   }
318  
319   delete fMatchedTrackIndex   ; 
320   delete fMatchedClusterIndex ; 
321   delete fResidualEta         ; 
322   delete fResidualPhi         ; 
323   delete fPIDUtils            ;
324
325   InitTrackCuts();
326 }
327
328 //_______________________________________________________________________________
329 Bool_t AliEMCALRecoUtils::AcceptCalibrateCell(const Int_t absID, const Int_t bc,
330                                               Float_t  & amp,    Double_t & time, 
331                                               AliVCaloCells* cells) 
332 {
333   // Reject cell if criteria not passed and calibrate it
334   
335   AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
336   
337   if(absID < 0 || absID >= 24*48*geom->GetNumberOfSuperModules()) return kFALSE;
338   
339   Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1; 
340   
341   if(!geom->GetCellIndex(absID,imod,iTower,iIphi,iIeta)) 
342   {
343     // cell absID does not exist
344     amp=0; time = 1.e9;
345     return kFALSE; 
346   }
347   
348   geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);        
349   
350   // Do not include bad channels found in analysis,
351   if( IsBadChannelsRemovalSwitchedOn() && GetEMCALChannelStatus(imod, ieta, iphi)) 
352   {
353     return kFALSE;
354   }
355   
356   //Recalibrate energy
357   amp  = cells->GetCellAmplitude(absID);
358   if(!fCellsRecalibrated && IsRecalibrationOn())
359     amp *= GetEMCALChannelRecalibrationFactor(imod,ieta,iphi);
360   
361   
362   // Recalibrate time
363   time = cells->GetCellTime(absID);
364   
365   RecalibrateCellTime(absID,bc,time);
366   
367   return kTRUE;
368 }
369
370 //_____________________________________________________________________________
371 Bool_t AliEMCALRecoUtils::CheckCellFiducialRegion(const AliEMCALGeometry* geom, 
372                                                   const AliVCluster* cluster, 
373                                                   AliVCaloCells* cells) 
374 {
375         // Given the list of AbsId of the cluster, get the maximum cell and 
376         // check if there are fNCellsFromBorder from the calorimeter border
377         
378   if(!cluster)
379   {
380     AliInfo("Cluster pointer null!");
381     return kFALSE;
382   }
383   
384   //If the distance to the border is 0 or negative just exit accept all clusters
385         if(cells->GetType()==AliVCaloCells::kEMCALCell && fNCellsFromEMCALBorder <= 0 ) return kTRUE;
386   
387   Int_t absIdMax        = -1, iSM =-1, ieta = -1, iphi = -1;
388   Bool_t shared = kFALSE;
389   GetMaxEnergyCell(geom, cells, cluster, absIdMax,  iSM, ieta, iphi, shared);
390   
391   AliDebug(2,Form("Cluster Max AbsId %d, Cell Energy %2.2f, Cluster Energy %2.2f, Ncells from border %d, EMCAL eta=0 %d\n", 
392                   absIdMax, cells->GetCellAmplitude(absIdMax), cluster->E(), fNCellsFromEMCALBorder, fNoEMCALBorderAtEta0));
393         
394         if(absIdMax==-1) return kFALSE;
395         
396         //Check if the cell is close to the borders:
397         Bool_t okrow = kFALSE;
398         Bool_t okcol = kFALSE;
399   
400   if(iSM < 0 || iphi < 0 || ieta < 0 ) 
401   {
402     AliFatal(Form("Negative value for super module: %d, or cell ieta: %d, or cell iphi: %d, check EMCAL geometry name\n",
403                   iSM,ieta,iphi));
404   }
405   
406   //Check rows/phi
407   if(iSM < 10)
408   {
409     if(iphi >= fNCellsFromEMCALBorder && iphi < 24-fNCellsFromEMCALBorder) okrow =kTRUE; 
410   }
411   else if (iSM >=10 && ( ( geom->GetEMCGeometry()->GetGeoName()).Contains("12SMV1"))) 
412   {
413     if(iphi >= fNCellsFromEMCALBorder && iphi < 8-fNCellsFromEMCALBorder) okrow =kTRUE; //1/3 sm case
414   }
415   else 
416   {
417     if(iphi >= fNCellsFromEMCALBorder && iphi < 12-fNCellsFromEMCALBorder) okrow =kTRUE; // half SM case
418   }
419   
420   //Check columns/eta
421   if(!fNoEMCALBorderAtEta0)
422   {
423     if(ieta  > fNCellsFromEMCALBorder && ieta < 48-fNCellsFromEMCALBorder) okcol =kTRUE; 
424   }
425   else
426   {
427     if(iSM%2==0)
428     {
429       if(ieta >= fNCellsFromEMCALBorder)     okcol = kTRUE;     
430     }
431     else 
432     {
433       if(ieta <  48-fNCellsFromEMCALBorder)  okcol = kTRUE;     
434     }
435   }//eta 0 not checked
436   
437   AliDebug(2,Form("EMCAL Cluster in %d cells fiducial volume: ieta %d, iphi %d, SM %d:  column? %d, row? %d\nq",
438                   fNCellsFromEMCALBorder, ieta, iphi, iSM, okcol, okrow));
439         
440         if (okcol && okrow) 
441   {
442     //printf("Accept\n");
443     return kTRUE;
444   }
445         else  
446   {
447     //printf("Reject\n");
448     AliDebug(2,Form("Reject cluster in border, max cell : ieta %d, iphi %d, SM %d\n",ieta, iphi, iSM));
449     return kFALSE;
450   }
451         
452 }       
453
454
455 //_______________________________________________________________________________
456 Bool_t AliEMCALRecoUtils::ClusterContainsBadChannel(const AliEMCALGeometry* geom, 
457                                                     const UShort_t* cellList, 
458                                                     const Int_t nCells)
459 {
460   // Check that in the cluster cells, there is no bad channel of those stored 
461   // in fEMCALBadChannelMap or fPHOSBadChannelMap
462         
463   if(!fRemoveBadChannels)  return kFALSE;
464   if(!fEMCALBadChannelMap) return kFALSE;
465         
466   Int_t icol = -1;
467   Int_t irow = -1;
468   Int_t imod = -1;
469   for(Int_t iCell = 0; iCell<nCells; iCell++)
470   {
471     //Get the column and row
472     Int_t iTower = -1, iIphi = -1, iIeta = -1; 
473     geom->GetCellIndex(cellList[iCell],imod,iTower,iIphi,iIeta); 
474     if(fEMCALBadChannelMap->GetEntries() <= imod) continue;
475     geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,irow,icol);                      
476     if(GetEMCALChannelStatus(imod, icol, irow))
477     {
478       AliDebug(2,Form("Cluster with bad channel: SM %d, col %d, row %d\n",imod, icol, irow));
479       return kTRUE;
480     }
481                 
482   }// cell cluster loop
483         
484   return kFALSE;
485 }
486
487 //_____________________________________________________________________________________________
488 Bool_t AliEMCALRecoUtils::IsExoticCell(const Int_t absID, AliVCaloCells* cells, const Int_t bc)
489 {
490   // Look to cell neighbourhood and reject if it seems exotic
491   // Do before recalibrating the cells
492
493   if(!fRejectExoticCells) return kFALSE;
494   
495   AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
496   
497   Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1; 
498   geom->GetCellIndex(absID,imod,iTower,iIphi,iIeta); 
499   geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);        
500   
501   //Get close cells index, energy and time, not in corners
502   
503   Int_t absID1 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi+1, ieta);
504   Int_t absID2 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi-1, ieta);
505   
506   // In case of cell in eta = 0 border, depending on SM shift the cross cell index
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, 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, iphi, AliEMCALGeoParams::fgkEMCALCols-1); 
519   }
520   else
521   {
522     absID3 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta+1);
523     absID4 = geom-> GetAbsCellIdFromCellIndexes(imod, iphi, ieta-1); 
524   }
525
526   Float_t  ecell  = 0, ecell1  = 0, ecell2  = 0, ecell3  = 0, ecell4  = 0;
527   Double_t tcell  = 0, tcell1  = 0, tcell2  = 0, tcell3  = 0, tcell4  = 0;
528   Bool_t   accept = 0, accept1 = 0, accept2 = 0, accept3 = 0, accept4 = 0;
529   
530   accept  = AcceptCalibrateCell(absID, bc, ecell ,tcell ,cells); 
531     
532   if(!accept) return kTRUE; // reject this cell
533   
534   if(ecell < fExoticCellMinAmplitude) return kFALSE; // do not reject low energy cells
535   
536   accept1 = AcceptCalibrateCell(absID1,bc, ecell1,tcell1,cells); 
537   accept2 = AcceptCalibrateCell(absID2,bc, ecell2,tcell2,cells); 
538   accept3 = AcceptCalibrateCell(absID3,bc, ecell3,tcell3,cells); 
539   accept4 = AcceptCalibrateCell(absID4,bc, ecell4,tcell4,cells); 
540   
541   /*
542     printf("Cell absID %d \n",absID);
543     printf("\t  accept1 %d, accept2 %d, accept3 %d, accept4 %d\n",
544            accept1,accept2,accept3,accept4);
545     printf("\t id %d: id1 %d, id2 %d, id3 %d, id4 %d\n",
546            absID,absID1,absID2,absID3,absID4);
547     printf("\t e %f: e1 %f, e2 %f, e3 %f, e4 %f\n",
548            ecell,ecell1,ecell2,ecell3,ecell4);
549     printf("\t t %f: t1 %f, t2 %f, t3 %f, t4 %f;\n dt1 %f, dt2 %f, dt3 %f, dt4 %f\n",
550            tcell*1.e9,tcell1*1.e9,tcell2*1.e9,tcell3*1.e9,tcell4*1.e9,
551            TMath::Abs(tcell-tcell1)*1.e9, TMath::Abs(tcell-tcell2)*1.e9, TMath::Abs(tcell-tcell3)*1.e9, TMath::Abs(tcell-tcell4)*1.e9);
552   */
553   
554   if(TMath::Abs(tcell-tcell1)*1.e9 > fExoticCellDiffTime) ecell1 = 0 ;
555   if(TMath::Abs(tcell-tcell2)*1.e9 > fExoticCellDiffTime) ecell2 = 0 ;
556   if(TMath::Abs(tcell-tcell3)*1.e9 > fExoticCellDiffTime) ecell3 = 0 ;
557   if(TMath::Abs(tcell-tcell4)*1.e9 > fExoticCellDiffTime) ecell4 = 0 ;
558
559   Float_t eCross = ecell1+ecell2+ecell3+ecell4;
560
561   //printf("\t eCell %f, eCross %f, 1-eCross/eCell %f\n",ecell,eCross,1-eCross/ecell);
562   
563   if(1-eCross/ecell > fExoticCellFraction) 
564   {
565     AliDebug(2,Form("AliEMCALRecoUtils::IsExoticCell() - EXOTIC CELL id %d, eCell %f, eCross %f, 1-eCross/eCell %f\n",
566                     absID,ecell,eCross,1-eCross/ecell));
567     return kTRUE;
568   }
569
570   return kFALSE;
571 }
572
573 //___________________________________________________________________
574 Bool_t AliEMCALRecoUtils::IsExoticCluster(const AliVCluster *cluster, 
575                                           AliVCaloCells *cells, 
576                                           const Int_t bc) 
577 {
578   // Check if the cluster highest energy tower is exotic
579   
580   if(!cluster)
581   {
582     AliInfo("Cluster pointer null!");
583     return kFALSE;
584   }
585   
586   if(!fRejectExoticCluster) return kFALSE;
587   
588   // Get highest energy tower
589   AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
590   Int_t iSupMod = -1, absId = -1, ieta = -1, iphi = -1;
591   Bool_t shared = kFALSE;
592   GetMaxEnergyCell(geom, cells, cluster, absId, iSupMod, ieta, iphi, shared);
593   
594   return IsExoticCell(absId,cells,bc);
595   
596 }
597
598 //_______________________________________________________________________
599 Float_t AliEMCALRecoUtils::SmearClusterEnergy(const AliVCluster* cluster) 
600 {
601   //In case of MC analysis, smear energy to match resolution/calibration in real data
602   
603   if(!cluster)
604   {
605     AliInfo("Cluster pointer null!");
606     return 0;
607   }
608   
609   Float_t energy    = cluster->E() ;
610   Float_t rdmEnergy = energy ;
611   if(fSmearClusterEnergy)
612   {
613     rdmEnergy = fRandom.Gaus(energy,fSmearClusterParam[0] * TMath::Sqrt(energy) +
614                                     fSmearClusterParam[1] * energy +
615                                     fSmearClusterParam[2] );
616     AliDebug(2, Form("Energy: original %f, smeared %f\n", energy, rdmEnergy));
617   }
618   
619   return rdmEnergy;
620 }
621
622 //____________________________________________________________________________
623 Float_t AliEMCALRecoUtils::CorrectClusterEnergyLinearity(AliVCluster* cluster)
624 {
625   // Correct cluster energy from non linearity functions
626   
627   if(!cluster)
628   {
629     AliInfo("Cluster pointer null!");
630     return 0;
631   }
632   
633   Float_t energy = cluster->E();
634
635   switch (fNonLinearityFunction) 
636   {
637       
638     case kPi0MC:
639     {
640       //Non-Linearity correction (from MC with function ([0]*exp(-[1]/E))+(([2]/([3]*2.*TMath::Pi())*exp(-(E-[4])^2/(2.*[3]^2)))))
641       //Double_t fNonLinearityParams[0] = 1.014;
642       //Double_t fNonLinearityParams[1] = -0.03329;
643       //Double_t fNonLinearityParams[2] = -0.3853;
644       //Double_t fNonLinearityParams[3] = 0.5423;
645       //Double_t fNonLinearityParams[4] = -0.4335;
646        energy *= (fNonLinearityParams[0]*exp(-fNonLinearityParams[1]/energy))+
647                   ((fNonLinearityParams[2]/(fNonLinearityParams[3]*2.*TMath::Pi())*
648                     exp(-(energy-fNonLinearityParams[4])*(energy-fNonLinearityParams[4])/(2.*fNonLinearityParams[3]*fNonLinearityParams[3]))));
649       break;
650     }
651      
652     case kPi0GammaGamma:
653     {
654       //Non-Linearity correction (from Olga Data with function p0+p1*exp(-p2*E))
655       //Double_t fNonLinearityParams[0] = 1.04;
656       //Double_t fNonLinearityParams[1] = -0.1445;
657       //Double_t fNonLinearityParams[2] = 1.046;
658       energy /= (fNonLinearityParams[0]+fNonLinearityParams[1]*exp(-fNonLinearityParams[2]*energy)); //Olga function
659       break;
660     }
661       
662     case kPi0GammaConversion:
663     {
664       //Non-Linearity correction (Nicolas from Dimitri Data with function C*[1-a*exp(-b*E)])
665       //fNonLinearityParams[0] = 0.139393/0.1349766;
666       //fNonLinearityParams[1] = 0.0566186;
667       //fNonLinearityParams[2] = 0.982133;
668       energy /= fNonLinearityParams[0]*(1-fNonLinearityParams[1]*exp(-fNonLinearityParams[2]*energy));
669       
670       break;
671     }
672       
673     case kBeamTest:
674     {
675       //From beam test, Alexei's results, for different ZS thresholds
676       //                        th=30 MeV; th = 45 MeV; th = 75 MeV
677       //fNonLinearityParams[0] = 1.007;      1.003;      1.002 
678       //fNonLinearityParams[1] = 0.894;      0.719;      0.797 
679       //fNonLinearityParams[2] = 0.246;      0.334;      0.358 
680       //Rescale the param[0] with 1.03
681       energy /= fNonLinearityParams[0]/(1+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]));
682       
683       break;
684     }
685       
686     case kBeamTestCorrected:
687     {
688       //From beam test, corrected for material between beam and EMCAL
689       //fNonLinearityParams[0] =  0.99078
690       //fNonLinearityParams[1] =  0.161499;
691       //fNonLinearityParams[2] =  0.655166; 
692       //fNonLinearityParams[3] =  0.134101;
693       //fNonLinearityParams[4] =  163.282;
694       //fNonLinearityParams[5] =  23.6904;
695       //fNonLinearityParams[6] =  0.978;
696         energy *= fNonLinearityParams[6]/(fNonLinearityParams[0]*(1./(1.+fNonLinearityParams[1]*exp(-energy/fNonLinearityParams[2]))*1./(1.+fNonLinearityParams[3]*exp((energy-fNonLinearityParams[4])/fNonLinearityParams[5]))));
697
698       break;
699     }
700       
701     case kNoCorrection:
702       AliDebug(2,"No correction on the energy\n");
703       break;
704       
705   }
706
707   return energy;
708 }
709
710 //__________________________________________________
711 void AliEMCALRecoUtils::InitNonLinearityParam()
712 {
713   //Initialising Non Linearity Parameters
714         
715   if(fNonLinearityFunction == kPi0MC) 
716   {
717     fNonLinearityParams[0] = 1.014;
718     fNonLinearityParams[1] = -0.03329;
719     fNonLinearityParams[2] = -0.3853;
720     fNonLinearityParams[3] = 0.5423;
721     fNonLinearityParams[4] = -0.4335;
722   }
723   
724   if(fNonLinearityFunction == kPi0GammaGamma) 
725   {
726     fNonLinearityParams[0] = 1.04;
727     fNonLinearityParams[1] = -0.1445;
728     fNonLinearityParams[2] = 1.046;
729   }     
730
731   if(fNonLinearityFunction == kPi0GammaConversion) 
732   {
733     fNonLinearityParams[0] = 0.139393;
734     fNonLinearityParams[1] = 0.0566186;
735     fNonLinearityParams[2] = 0.982133;
736   }     
737
738   if(fNonLinearityFunction == kBeamTest) 
739   {
740     if(fNonLinearThreshold == 30) 
741     {
742       fNonLinearityParams[0] = 1.007; 
743       fNonLinearityParams[1] = 0.894; 
744       fNonLinearityParams[2] = 0.246; 
745     }
746     if(fNonLinearThreshold == 45) 
747     {
748       fNonLinearityParams[0] = 1.003; 
749       fNonLinearityParams[1] = 0.719; 
750       fNonLinearityParams[2] = 0.334; 
751     }
752     if(fNonLinearThreshold == 75) 
753     {
754       fNonLinearityParams[0] = 1.002; 
755       fNonLinearityParams[1] = 0.797; 
756       fNonLinearityParams[2] = 0.358; 
757     }
758   }
759
760   if(fNonLinearityFunction == kBeamTestCorrected) 
761   {
762     fNonLinearityParams[0] =  0.99078;
763     fNonLinearityParams[1] =  0.161499;
764     fNonLinearityParams[2] =  0.655166; 
765     fNonLinearityParams[3] =  0.134101;
766     fNonLinearityParams[4] =  163.282;
767     fNonLinearityParams[5] =  23.6904;
768     fNonLinearityParams[6] =  0.978;
769   }
770 }
771
772 //_________________________________________________________
773 Float_t  AliEMCALRecoUtils::GetDepth(const Float_t energy, 
774                                      const Int_t iParticle, 
775                                      const Int_t iSM) const 
776 {
777   //Calculate shower depth for a given cluster energy and particle type
778
779   // parameters 
780   Float_t x0    = 1.31;
781   Float_t ecr   = 8;
782   Float_t depth = 0;
783   
784   switch ( iParticle )
785   {
786     case kPhoton:
787       depth = x0 * (TMath::Log(energy*1000/ ecr) + 0.5); //Multiply energy by 1000 to transform to MeV
788       break;
789       
790     case kElectron:
791       depth = x0 * (TMath::Log(energy*1000/ ecr) - 0.5); //Multiply energy by 1000 to transform to MeV
792       break;
793       
794     case kHadron:
795       // hadron 
796       // boxes anc. here
797       if(gGeoManager)
798       {
799         gGeoManager->cd("ALIC_1/XEN1_1");
800         TGeoNode        *geoXEn1    = gGeoManager->GetCurrentNode();
801         TGeoNodeMatrix  *geoSM      = dynamic_cast<TGeoNodeMatrix *>(geoXEn1->GetDaughter(iSM));
802         if(geoSM)
803         {
804           TGeoVolume      *geoSMVol   = geoSM->GetVolume(); 
805           TGeoShape       *geoSMShape = geoSMVol->GetShape();
806           TGeoBBox        *geoBox     = dynamic_cast<TGeoBBox *>(geoSMShape);
807           if(geoBox) depth = 0.5 * geoBox->GetDX()*2 ;
808           else AliFatal("Null GEANT box");
809         }
810         else AliFatal("NULL  GEANT node matrix");
811       }
812       else
813       {//electron
814         depth = x0 * (TMath::Log(energy*1000 / ecr)  - 0.5); //Multiply energy by 1000 to transform to MeV
815       }
816         
817       break;
818       
819     default://photon
820       depth = x0 * (TMath::Log(energy*1000 / ecr) + 0.5); //Multiply energy by 1000 to transform to MeV
821   }  
822   
823   return depth;
824 }
825
826 //____________________________________________________________________
827 void AliEMCALRecoUtils::GetMaxEnergyCell(const AliEMCALGeometry *geom, 
828                                          AliVCaloCells* cells, 
829                                          const AliVCluster* clu, 
830                                          Int_t  & absId,  
831                                          Int_t  & iSupMod, 
832                                          Int_t  & ieta, 
833                                          Int_t  & iphi, 
834                                          Bool_t & shared)
835 {
836   //For a given CaloCluster gets the absId of the cell 
837   //with maximum energy deposit.
838   
839   Double_t eMax        = -1.;
840   Double_t eCell       = -1.;
841   Float_t  fraction    = 1.;
842   Float_t  recalFactor = 1.;
843   Int_t    cellAbsId   = -1 ;
844
845   Int_t iTower  = -1;
846   Int_t iIphi   = -1;
847   Int_t iIeta   = -1;
848   Int_t iSupMod0= -1;
849
850   if(!clu)
851   {
852     AliInfo("Cluster pointer null!");
853     absId=-1; iSupMod0=-1, ieta = -1; iphi = -1; shared = -1;
854     return;
855   }
856   
857   for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++) 
858   {
859     cellAbsId = clu->GetCellAbsId(iDig);
860     fraction  = clu->GetCellAmplitudeFraction(iDig);
861     //printf("a Cell %d, id, %d, amp %f, fraction %f\n",iDig,cellAbsId,cells->GetCellAmplitude(cellAbsId),fraction);
862     if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
863     geom->GetCellIndex(cellAbsId,iSupMod,iTower,iIphi,iIeta); 
864     geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);
865     if     (iDig==0) 
866     {
867       iSupMod0=iSupMod;
868     }
869     else if(iSupMod0!=iSupMod) 
870     {
871       shared = kTRUE;
872       //printf("AliEMCALRecoUtils::GetMaxEnergyCell() - SHARED CLUSTER\n");
873     }
874     if(!fCellsRecalibrated && IsRecalibrationOn()) 
875     {
876       recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
877     }
878     eCell  = cells->GetCellAmplitude(cellAbsId)*fraction*recalFactor;
879     //printf("b Cell %d, id, %d, amp %f, fraction %f\n",iDig,cellAbsId,eCell,fraction);
880     if(eCell > eMax)  
881     { 
882       eMax  = eCell; 
883       absId = cellAbsId;
884       //printf("\t new max: cell %d, e %f, ecell %f\n",maxId, eMax,eCell);
885     }
886   }// cell loop
887   
888   //Get from the absid the supermodule, tower and eta/phi numbers
889   geom->GetCellIndex(absId,iSupMod,iTower,iIphi,iIeta); 
890   //Gives SuperModule and Tower numbers
891   geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,
892                                          iIphi, iIeta,iphi,ieta); 
893   //printf("Max id %d, iSM %d, col %d, row %d\n",absId,iSupMod,ieta,iphi);
894   //printf("Max end---\n");
895 }
896
897 //______________________________________
898 void AliEMCALRecoUtils::InitParameters()
899 {
900   // Initialize data members with default values
901   
902   fParticleType = kPhoton;
903   fPosAlgo      = kUnchanged;
904   fW0           = 4.5;
905   
906   fNonLinearityFunction = kNoCorrection;
907   fNonLinearThreshold   = 30;
908   
909   fExoticCellFraction     = 0.97;
910   fExoticCellDiffTime     = 1e6;
911   fExoticCellMinAmplitude = 0.5;
912   
913   fAODFilterMask = 32;
914   
915   fCutEtaPhiSum      = kTRUE;
916   fCutEtaPhiSeparate = kFALSE;
917   
918   fCutR   = 0.05; 
919   fCutEta = 0.025; 
920   fCutPhi = 0.05;
921   
922   fClusterWindow = 100;
923   fMass          = 0.139;
924   
925   fStepSurface   = 20.;                      
926   fStepCluster   = 5.;
927   fTrackCutsType = kLooseCut;
928   
929   fCutMinTrackPt     = 0;
930   fCutMinNClusterTPC = -1;
931   fCutMinNClusterITS = -1;
932   
933   fCutMaxChi2PerClusterTPC  = 1e10;
934   fCutMaxChi2PerClusterITS  = 1e10;
935   
936   fCutRequireTPCRefit     = kFALSE;
937   fCutRequireITSRefit     = kFALSE;
938   fCutAcceptKinkDaughters = kFALSE;
939   
940   fCutMaxDCAToVertexXY = 1e10;             
941   fCutMaxDCAToVertexZ  = 1e10;              
942   fCutDCAToVertex2D    = kFALSE;
943   
944   
945   //Misalignment matrices
946   for(Int_t i = 0; i < 15 ; i++) 
947   {
948     fMisalTransShift[i] = 0.; 
949     fMisalRotShift[i]   = 0.; 
950   }
951   
952   //Non linearity
953   for(Int_t i = 0; i < 7  ; i++) fNonLinearityParams[i] = 0.; 
954   
955   //For kBeamTestCorrected case, but default is no correction
956   fNonLinearityParams[0] =  0.99078;
957   fNonLinearityParams[1] =  0.161499;
958   fNonLinearityParams[2] =  0.655166; 
959   fNonLinearityParams[3] =  0.134101;
960   fNonLinearityParams[4] =  163.282;
961   fNonLinearityParams[5] =  23.6904;
962   fNonLinearityParams[6] =  0.978;
963   
964   //For kPi0GammaGamma case
965   //fNonLinearityParams[0] = 0.1457/0.1349766/1.038;
966   //fNonLinearityParams[1] = -0.02024/0.1349766/1.038;
967   //fNonLinearityParams[2] = 1.046;
968   
969   //Cluster energy smearing
970   fSmearClusterEnergy   = kFALSE;
971   fSmearClusterParam[0] = 0.07; // * sqrt E term
972   fSmearClusterParam[1] = 0.00; // * E term
973   fSmearClusterParam[2] = 0.00; // constant
974 }
975
976 //_____________________________________________________
977 void AliEMCALRecoUtils::InitEMCALRecalibrationFactors()
978 {
979   //Init EMCAL recalibration factors
980   AliDebug(2,"AliCalorimeterUtils::InitEMCALRecalibrationFactors()");
981         //In order to avoid rewriting the same histograms
982   Bool_t oldStatus = TH1::AddDirectoryStatus();
983   TH1::AddDirectory(kFALSE);
984   
985   fEMCALRecalibrationFactors = new TObjArray(12);
986   for (int i = 0; i < 12; i++) 
987     fEMCALRecalibrationFactors->Add(new TH2F(Form("EMCALRecalFactors_SM%d",i),
988                                              Form("EMCALRecalFactors_SM%d",i),  48, 0, 48, 24, 0, 24));
989   //Init the histograms with 1
990   for (Int_t sm = 0; sm < 12; sm++) 
991   {
992     for (Int_t i = 0; i < 48; i++) 
993     {
994       for (Int_t j = 0; j < 24; j++) 
995       {
996         SetEMCALChannelRecalibrationFactor(sm,i,j,1.);
997       }
998     }
999   }
1000   
1001   fEMCALRecalibrationFactors->SetOwner(kTRUE);
1002   fEMCALRecalibrationFactors->Compress();
1003         
1004   //In order to avoid rewriting the same histograms
1005   TH1::AddDirectory(oldStatus);         
1006 }
1007
1008 //_________________________________________________________
1009 void AliEMCALRecoUtils::InitEMCALTimeRecalibrationFactors()
1010 {
1011   //Init EMCAL recalibration factors
1012   AliDebug(2,"AliCalorimeterUtils::InitEMCALRecalibrationFactors()");
1013   //In order to avoid rewriting the same histograms
1014   Bool_t oldStatus = TH1::AddDirectoryStatus();
1015   TH1::AddDirectory(kFALSE);
1016   
1017   fEMCALTimeRecalibrationFactors = new TObjArray(4);
1018   for (int i = 0; i < 4; i++) 
1019     fEMCALTimeRecalibrationFactors->Add(new TH1F(Form("hAllTimeAvBC%d",i),
1020                                                  Form("hAllTimeAvBC%d",i),  
1021                                                  48*24*12,0.,48*24*12)          );
1022   //Init the histograms with 1
1023   for (Int_t bc = 0; bc < 4; bc++) 
1024   {
1025     for (Int_t i = 0; i < 48*24*12; i++) 
1026       SetEMCALChannelTimeRecalibrationFactor(bc,i,0.);
1027   }
1028   
1029   fEMCALTimeRecalibrationFactors->SetOwner(kTRUE);
1030   fEMCALTimeRecalibrationFactors->Compress();
1031   
1032   //In order to avoid rewriting the same histograms
1033   TH1::AddDirectory(oldStatus);         
1034 }
1035
1036 //____________________________________________________
1037 void AliEMCALRecoUtils::InitEMCALBadChannelStatusMap()
1038 {
1039   //Init EMCAL bad channels map
1040   AliDebug(2,"AliEMCALRecoUtils::InitEMCALBadChannelStatusMap()");
1041   //In order to avoid rewriting the same histograms
1042   Bool_t oldStatus = TH1::AddDirectoryStatus();
1043   TH1::AddDirectory(kFALSE);
1044   
1045   fEMCALBadChannelMap = new TObjArray(12);
1046   //TH2F * hTemp = new  TH2I("EMCALBadChannelMap","EMCAL SuperModule bad channel map", 48, 0, 48, 24, 0, 24);
1047   for (int i = 0; i < 12; i++) 
1048   {
1049     fEMCALBadChannelMap->Add(new TH2I(Form("EMCALBadChannelMap_Mod%d",i),Form("EMCALBadChannelMap_Mod%d",i), 48, 0, 48, 24, 0, 24));
1050   }
1051   
1052   fEMCALBadChannelMap->SetOwner(kTRUE);
1053   fEMCALBadChannelMap->Compress();
1054   
1055   //In order to avoid rewriting the same histograms
1056   TH1::AddDirectory(oldStatus);         
1057 }
1058
1059 //____________________________________________________________________________
1060 void AliEMCALRecoUtils::RecalibrateClusterEnergy(const AliEMCALGeometry* geom, 
1061                                                  AliVCluster * cluster, 
1062                                                  AliVCaloCells * cells, 
1063                                                  const Int_t bc)
1064 {
1065   // Recalibrate the cluster energy and Time, considering the recalibration map 
1066   // and the energy of the cells and time that compose the cluster.
1067   // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1068         
1069   if(!cluster)
1070   {
1071     AliInfo("Cluster pointer null!");
1072     return;
1073   }  
1074   
1075   //Get the cluster number of cells and list of absId, check what kind of cluster do we have.
1076   UShort_t * index    = cluster->GetCellsAbsId() ;
1077   Double_t * fraction = cluster->GetCellsAmplitudeFraction() ;
1078   Int_t ncells = cluster->GetNCells();
1079         
1080   //Initialize some used variables
1081   Float_t energy = 0;
1082   Int_t   absId  =-1;
1083   Int_t   icol   =-1, irow =-1, imod=1;
1084   Float_t factor = 1, frac = 0;
1085   Int_t   absIdMax = -1;
1086   Float_t emax     = 0;
1087   
1088   //Loop on the cells, get the cell amplitude and recalibration factor, multiply and and to the new energy
1089   for(Int_t icell = 0; icell < ncells; icell++)
1090   {
1091     absId = index[icell];
1092     frac =  fraction[icell];
1093     if(frac < 1e-5) frac = 1; //in case of EMCAL, this is set as 0 since unfolding is off
1094     
1095     if(!fCellsRecalibrated && IsRecalibrationOn()) 
1096     {
1097       // Energy  
1098       Int_t iTower = -1, iIphi = -1, iIeta = -1; 
1099       geom->GetCellIndex(absId,imod,iTower,iIphi,iIeta); 
1100       if(fEMCALRecalibrationFactors->GetEntries() <= imod) continue;
1101       geom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,irow,icol);                    
1102       factor = GetEMCALChannelRecalibrationFactor(imod,icol,irow);
1103       
1104       AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - recalibrate cell: module %d, col %d, row %d, cell fraction %f,recalibration factor %f, cell energy %f\n",
1105                       imod,icol,irow,frac,factor,cells->GetCellAmplitude(absId)));
1106       
1107     } 
1108     
1109     energy += cells->GetCellAmplitude(absId)*factor*frac;
1110     
1111     if(emax < cells->GetCellAmplitude(absId)*factor*frac)
1112     {
1113       emax     = cells->GetCellAmplitude(absId)*factor*frac;
1114       absIdMax = absId;
1115     }
1116   }
1117         
1118   AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - Energy before %f, after %f \n",cluster->E(),energy));
1119
1120   cluster->SetE(energy);
1121
1122   // Recalculate time of cluster   
1123   Double_t timeorg = cluster->GetTOF();
1124   if(!fCellsRecalibrated && IsTimeRecalibrationOn())
1125   {
1126     Double_t time = timeorg;
1127     RecalibrateCellTime(absIdMax,bc,time);
1128     cluster->SetTOF(time);
1129   } 
1130
1131   AliDebug(2,Form("AliEMCALRecoUtils::RecalibrateClusterEnergy - Time before %f, after %f \n",timeorg,cluster->GetTOF()));
1132
1133 }
1134
1135 //_____________________________________________________________
1136 void AliEMCALRecoUtils::RecalibrateCells(AliVCaloCells * cells,
1137                                          const Int_t bc)
1138 {
1139   // Recalibrate the cells time and energy, considering the recalibration map and the energy 
1140   // of the cells that compose the cluster.
1141   // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1142
1143   if(!IsRecalibrationOn() && !IsTimeRecalibrationOn()) return;
1144   
1145   if(!cells)
1146   {
1147     AliInfo("Cells pointer null!");
1148     return;
1149   }  
1150   
1151   Int_t    absId  =-1;
1152   Bool_t   accept = kFALSE;
1153   Float_t  ecell  = 0;
1154   Double_t tcell  = 0;
1155   
1156   Int_t nEMcell  = cells->GetNumberOfCells() ;  
1157   for (Int_t iCell = 0; iCell < nEMcell; iCell++) 
1158   { 
1159     absId  = cells->GetCellNumber(iCell);
1160     
1161     accept = AcceptCalibrateCell(absId, bc, ecell ,tcell ,cells); 
1162     if(!accept) 
1163     {
1164       ecell = 0;
1165       tcell = 0;
1166     }
1167     
1168     //Set new values
1169     cells->SetCell(iCell,absId,ecell, tcell);
1170   }
1171
1172   fCellsRecalibrated = kTRUE;
1173 }
1174
1175 //_______________________________________________________________________________________________________
1176 void AliEMCALRecoUtils::RecalibrateCellTime(const Int_t absId, const Int_t bc, Double_t & celltime) const
1177 {
1178   // Recalibrate time of cell with absID  considering the recalibration map 
1179   // bc= bunch crossing number returned by esdevent->GetBunchCrossNumber();
1180   
1181   if(!fCellsRecalibrated && IsTimeRecalibrationOn() && bc >= 0)
1182   {
1183     celltime -= GetEMCALChannelTimeRecalibrationFactor(bc%4,absId)*1.e-9;    ;  
1184   }
1185 }
1186   
1187 //______________________________________________________________________________
1188 void AliEMCALRecoUtils::RecalculateClusterPosition(const AliEMCALGeometry *geom, 
1189                                                    AliVCaloCells* cells, 
1190                                                    AliVCluster* clu)
1191 {
1192   //For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1193   
1194   if(!clu)
1195   {
1196     AliInfo("Cluster pointer null!");
1197     return;
1198   }
1199     
1200   if     (fPosAlgo==kPosTowerGlobal) RecalculateClusterPositionFromTowerGlobal( geom, cells, clu);
1201   else if(fPosAlgo==kPosTowerIndex)  RecalculateClusterPositionFromTowerIndex ( geom, cells, clu);
1202   else   AliDebug(2,"Algorithm to recalculate position not selected, do nothing.");
1203 }  
1204
1205 //_____________________________________________________________________________________________
1206 void AliEMCALRecoUtils::RecalculateClusterPositionFromTowerGlobal(const AliEMCALGeometry *geom, 
1207                                                                   AliVCaloCells* cells, 
1208                                                                   AliVCluster* clu)
1209 {
1210   // For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1211   // The algorithm is a copy of what is done in AliEMCALRecPoint
1212   
1213   Double_t eCell       = 0.;
1214   Float_t  fraction    = 1.;
1215   Float_t  recalFactor = 1.;
1216   
1217   Int_t    absId   = -1;
1218   Int_t    iTower  = -1, iIphi  = -1, iIeta  = -1;
1219   Int_t    iSupModMax = -1, iSM=-1, iphi   = -1, ieta   = -1;
1220   Float_t  weight = 0.,  totalWeight=0.;
1221   Float_t  newPos[3] = {0,0,0};
1222   Double_t pLocal[3], pGlobal[3];
1223   Bool_t shared = kFALSE;
1224
1225   Float_t  clEnergy = clu->E(); //Energy already recalibrated previously
1226   GetMaxEnergyCell(geom, cells, clu, absId,  iSupModMax, ieta, iphi,shared);
1227   Double_t depth = GetDepth(clEnergy,fParticleType,iSupModMax) ;
1228   
1229   //printf("** Cluster energy %f, ncells %d, depth %f\n",clEnergy,clu->GetNCells(),depth);
1230   
1231   for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++) 
1232   {
1233     absId = clu->GetCellAbsId(iDig);
1234     fraction  = clu->GetCellAmplitudeFraction(iDig);
1235     if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1236     
1237     if (!fCellsRecalibrated)
1238     {
1239       geom->GetCellIndex(absId,iSM,iTower,iIphi,iIeta); 
1240       geom->GetCellPhiEtaIndexInSModule(iSM,iTower,iIphi, iIeta,iphi,ieta);                     
1241       
1242       if(IsRecalibrationOn()) 
1243       {
1244         recalFactor = GetEMCALChannelRecalibrationFactor(iSM,ieta,iphi);
1245       }
1246     }
1247     
1248     eCell  = cells->GetCellAmplitude(absId)*fraction*recalFactor;
1249     
1250     weight = GetCellWeight(eCell,clEnergy);
1251     totalWeight += weight;
1252     
1253     geom->RelPosCellInSModule(absId,depth,pLocal[0],pLocal[1],pLocal[2]);
1254     //printf("pLocal (%f,%f,%f), SM %d, absId %d\n",pLocal[0],pLocal[1],pLocal[2],iSupModMax,absId);
1255     geom->GetGlobal(pLocal,pGlobal,iSupModMax);
1256     //printf("pLocal (%f,%f,%f)\n",pGlobal[0],pGlobal[1],pGlobal[2]);
1257
1258     for(int i=0; i<3; i++ ) newPos[i] += (weight*pGlobal[i]);
1259   }// cell loop
1260   
1261   if(totalWeight>0)
1262   {
1263     for(int i=0; i<3; i++ )    newPos[i] /= totalWeight;
1264   }
1265     
1266   //Float_t pos[]={0,0,0};
1267   //clu->GetPosition(pos);
1268   //printf("OldPos  : %2.3f,%2.3f,%2.3f\n",pos[0],pos[1],pos[2]);
1269   //printf("NewPos  : %2.3f,%2.3f,%2.3f\n",newPos[0],newPos[1],newPos[2]);
1270   
1271   if(iSupModMax > 1) //sector 1
1272   {
1273     newPos[0] +=fMisalTransShift[3];//-=3.093; 
1274     newPos[1] +=fMisalTransShift[4];//+=6.82;
1275     newPos[2] +=fMisalTransShift[5];//+=1.635;
1276     //printf("   +    : %2.3f,%2.3f,%2.3f\n",fMisalTransShift[3],fMisalTransShift[4],fMisalTransShift[5]);
1277   } else //sector 0
1278   {
1279     newPos[0] +=fMisalTransShift[0];//+=1.134;
1280     newPos[1] +=fMisalTransShift[1];//+=8.2;
1281     newPos[2] +=fMisalTransShift[2];//+=1.197;
1282     //printf("   +    : %2.3f,%2.3f,%2.3f\n",fMisalTransShift[0],fMisalTransShift[1],fMisalTransShift[2]);
1283   }
1284   //printf("NewPos : %2.3f,%2.3f,%2.3f\n",newPos[0],newPos[1],newPos[2]);
1285
1286   clu->SetPosition(newPos);
1287 }  
1288
1289 //____________________________________________________________________________________________
1290 void AliEMCALRecoUtils::RecalculateClusterPositionFromTowerIndex(const AliEMCALGeometry *geom, 
1291                                                                  AliVCaloCells* cells, 
1292                                                                  AliVCluster* clu)
1293 {
1294   // For a given CaloCluster recalculates the position for a given set of misalignment shifts and puts it again in the CaloCluster.
1295   // The algorithm works with the tower indeces, averages the indeces and from them it calculates the global position
1296   
1297   Double_t eCell       = 1.;
1298   Float_t  fraction    = 1.;
1299   Float_t  recalFactor = 1.;
1300   
1301   Int_t absId   = -1;
1302   Int_t iTower  = -1;
1303   Int_t iIphi   = -1, iIeta   = -1;
1304   Int_t iSupMod = -1, iSupModMax = -1;
1305   Int_t iphi = -1, ieta =-1;
1306   Bool_t shared = kFALSE;
1307
1308   Float_t clEnergy = clu->E(); //Energy already recalibrated previously.
1309   GetMaxEnergyCell(geom, cells, clu, absId,  iSupModMax, ieta, iphi,shared);
1310   Float_t  depth = GetDepth(clEnergy,fParticleType,iSupMod) ;
1311
1312   Float_t weight = 0., weightedCol = 0., weightedRow = 0., totalWeight=0.;
1313   Bool_t areInSameSM = kTRUE; //exclude clusters with cells in different SMs for now
1314   Int_t startingSM = -1;
1315   
1316   for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++) 
1317   {
1318     absId = clu->GetCellAbsId(iDig);
1319     fraction  = clu->GetCellAmplitudeFraction(iDig);
1320     if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1321
1322     if     (iDig==0)  startingSM = iSupMod;
1323     else if(iSupMod != startingSM) areInSameSM = kFALSE;
1324
1325     eCell  = cells->GetCellAmplitude(absId);
1326     
1327     geom->GetCellIndex(absId,iSupMod,iTower,iIphi,iIeta); 
1328     geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);           
1329     
1330     if (!fCellsRecalibrated)
1331     {
1332       if(IsRecalibrationOn()) 
1333       {
1334         recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1335       }
1336     }
1337     
1338     eCell  = cells->GetCellAmplitude(absId)*fraction*recalFactor;
1339     
1340     weight = GetCellWeight(eCell,clEnergy);
1341     if(weight < 0) weight = 0;
1342     totalWeight += weight;
1343     weightedCol += ieta*weight;
1344     weightedRow += iphi*weight;
1345     
1346     //printf("Max cell? cell %d, amplitude org %f, fraction %f, recalibration %f, amplitude new %f \n",cellAbsId, cells->GetCellAmplitude(cellAbsId), fraction, recalFactor, eCell) ;
1347   }// cell loop
1348     
1349   Float_t xyzNew[]={0.,0.,0.};
1350   if(areInSameSM == kTRUE) 
1351   {
1352     //printf("In Same SM\n");
1353     weightedCol = weightedCol/totalWeight;
1354     weightedRow = weightedRow/totalWeight;
1355     geom->RecalculateTowerPosition(weightedRow, weightedCol, iSupModMax, depth, fMisalTransShift, fMisalRotShift, xyzNew); 
1356   } 
1357   else 
1358   {
1359     //printf("In Different SM\n");
1360     geom->RecalculateTowerPosition(iphi,        ieta,        iSupModMax, depth, fMisalTransShift, fMisalRotShift, xyzNew); 
1361   }
1362   
1363   clu->SetPosition(xyzNew);
1364 }
1365
1366 //___________________________________________________________________________________________
1367 void AliEMCALRecoUtils::RecalculateClusterDistanceToBadChannel(const AliEMCALGeometry * geom, 
1368                                                                AliVCaloCells* cells, 
1369                                                                AliVCluster * cluster)
1370 {           
1371   //re-evaluate distance to bad channel with updated bad map
1372   
1373   if(!fRecalDistToBadChannels) return;
1374   
1375   if(!cluster)
1376   {
1377     AliInfo("Cluster pointer null!");
1378     return;
1379   }  
1380   
1381   //Get channels map of the supermodule where the cluster is.
1382   Int_t absIdMax        = -1, iSupMod =-1, icolM = -1, irowM = -1;
1383   Bool_t shared = kFALSE;
1384   GetMaxEnergyCell(geom, cells, cluster, absIdMax,  iSupMod, icolM, irowM, shared);
1385   TH2D* hMap  = (TH2D*)fEMCALBadChannelMap->At(iSupMod);
1386
1387   Int_t dRrow, dRcol;   
1388   Float_t  minDist = 10000.;
1389   Float_t  dist    = 0.;
1390   
1391   //Loop on tower status map 
1392   for(Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++)
1393   {
1394     for(Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++)
1395     {
1396       //Check if tower is bad.
1397       if(hMap->GetBinContent(icol,irow)==0) continue;
1398       //printf("AliEMCALRecoUtils::RecalculateDistanceToBadChannels() - \n \t Bad channel in SM %d, col %d, row %d, \n \t Cluster max in col %d, row %d\n",
1399       //       iSupMod,icol, irow, icolM,irowM);
1400       
1401       dRrow=TMath::Abs(irowM-irow);
1402       dRcol=TMath::Abs(icolM-icol);
1403       dist=TMath::Sqrt(dRrow*dRrow+dRcol*dRcol);
1404       if(dist < minDist)
1405       {
1406         //printf("MIN DISTANCE TO BAD %2.2f\n",dist);
1407         minDist = dist;
1408       }
1409     }
1410   }
1411   
1412   //In case the cluster is shared by 2 SuperModules, need to check the map of the second Super Module
1413   if (shared) 
1414   {
1415     TH2D* hMap2 = 0;
1416     Int_t iSupMod2 = -1;
1417     
1418     //The only possible combinations are (0,1), (2,3) ... (8,9)
1419     if(iSupMod%2) iSupMod2 = iSupMod-1;
1420     else          iSupMod2 = iSupMod+1;
1421     hMap2  = (TH2D*)fEMCALBadChannelMap->At(iSupMod2);
1422     
1423     //Loop on tower status map of second super module
1424     for(Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++)
1425     {
1426       for(Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++)
1427       {
1428         //Check if tower is bad.
1429         if(hMap2->GetBinContent(icol,irow)==0) continue;
1430         //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",
1431         //     iSupMod2,icol, irow,iSupMod,icolM,irowM);
1432         dRrow=TMath::Abs(irow-irowM);
1433         
1434         if(iSupMod%2) 
1435         {
1436           dRcol=TMath::Abs(icol-(AliEMCALGeoParams::fgkEMCALCols+icolM));
1437         } else 
1438         {
1439           dRcol=TMath::Abs(AliEMCALGeoParams::fgkEMCALCols+icol-icolM);
1440         }                    
1441         
1442         dist=TMath::Sqrt(dRrow*dRrow+dRcol*dRcol);
1443         if(dist < minDist) minDist = dist;        
1444       }
1445     }
1446   }// shared cluster in 2 SuperModules
1447   
1448   AliDebug(2,Form("Max cluster cell (SM,col,row)=(%d %d %d) - Distance to Bad Channel %2.2f",iSupMod, icolM, irowM, minDist));
1449   cluster->SetDistanceToBadChannel(minDist);
1450 }
1451
1452 //__________________________________________________________________
1453 void AliEMCALRecoUtils::RecalculateClusterPID(AliVCluster * cluster)
1454 {           
1455   //re-evaluate identification parameters with bayesian
1456   
1457   if(!cluster)
1458   {
1459     AliInfo("Cluster pointer null!");
1460     return;
1461   }
1462   
1463   if ( cluster->GetM02() != 0)
1464     fPIDUtils->ComputePID(cluster->E(),cluster->GetM02());
1465   
1466   Float_t pidlist[AliPID::kSPECIESN+1];
1467         for(Int_t i = 0; i < AliPID::kSPECIESN+1; i++) pidlist[i] = fPIDUtils->GetPIDFinal(i);
1468         
1469   cluster->SetPID(pidlist);
1470 }
1471
1472 //___________________________________________________________________________________________________________________
1473 void AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(const AliEMCALGeometry * geom, 
1474                                                                 AliVCaloCells* cells, 
1475                                                                 AliVCluster * cluster,
1476                                                                 Float_t & l0,   Float_t & l1,   
1477                                                                 Float_t & disp, Float_t & dEta, Float_t & dPhi,
1478                                                                 Float_t & sEta, Float_t & sPhi, Float_t & sEtaPhi)
1479 {
1480   // Calculates new center of gravity in the local EMCAL-module coordinates 
1481   // and tranfers into global ALICE coordinates
1482   // Calculates Dispersion and main axis
1483   
1484   if(!cluster)
1485   {
1486     AliInfo("Cluster pointer null!");
1487     return;
1488   }
1489     
1490   Double_t eCell       = 0.;
1491   Float_t  fraction    = 1.;
1492   Float_t  recalFactor = 1.;
1493
1494   Int_t    iSupMod = -1;
1495   Int_t    iTower  = -1;
1496   Int_t    iIphi   = -1;
1497   Int_t    iIeta   = -1;
1498   Int_t    iphi    = -1;
1499   Int_t    ieta    = -1;
1500   Double_t etai    = -1.;
1501   Double_t phii    = -1.;
1502   
1503   Int_t    nstat   = 0 ;
1504   Float_t  wtot    = 0.;
1505   Double_t w       = 0.;
1506   Double_t etaMean = 0.;
1507   Double_t phiMean = 0.;
1508     
1509   //Loop on cells
1510   for(Int_t iDigit=0; iDigit < cluster->GetNCells(); iDigit++) 
1511   {
1512     //Get from the absid the supermodule, tower and eta/phi numbers
1513     geom->GetCellIndex(cluster->GetCellAbsId(iDigit),iSupMod,iTower,iIphi,iIeta); 
1514     geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi,iIeta, iphi,ieta);        
1515     
1516     //Get the cell energy, if recalibration is on, apply factors
1517     fraction  = cluster->GetCellAmplitudeFraction(iDigit);
1518     if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1519     
1520     if (!fCellsRecalibrated)
1521     {
1522       if(IsRecalibrationOn()) 
1523       {
1524         recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1525       }
1526     }
1527     
1528     eCell  = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1529     
1530     if(cluster->E() > 0 && eCell > 0)
1531     {
1532       w  = GetCellWeight(eCell,cluster->E());
1533       
1534       etai=(Double_t)ieta;
1535       phii=(Double_t)iphi;      
1536       
1537       if(w > 0.0) 
1538       {
1539         wtot += w ;
1540         nstat++;                        
1541         //Shower shape
1542         sEta     += w * etai * etai ;
1543         etaMean  += w * etai ;
1544         sPhi     += w * phii * phii ;
1545         phiMean  += w * phii ; 
1546         sEtaPhi  += w * etai * phii ; 
1547       }
1548     } 
1549     else
1550       AliError(Form("Wrong energy %f and/or amplitude %f\n", eCell, cluster->E()));
1551   }//cell loop
1552   
1553   //Normalize to the weight     
1554   if (wtot > 0) 
1555   {
1556     etaMean /= wtot ;
1557     phiMean /= wtot ;
1558   }
1559   else
1560     AliError(Form("Wrong weight %f\n", wtot));
1561   
1562   //Calculate dispersion        
1563   for(Int_t iDigit=0; iDigit < cluster->GetNCells(); iDigit++) 
1564   {
1565     //Get from the absid the supermodule, tower and eta/phi numbers
1566     geom->GetCellIndex(cluster->GetCellAbsId(iDigit),iSupMod,iTower,iIphi,iIeta); 
1567     geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi,iIeta, iphi,ieta);
1568     
1569     //Get the cell energy, if recalibration is on, apply factors
1570     fraction  = cluster->GetCellAmplitudeFraction(iDigit);
1571     if(fraction < 1e-4) fraction = 1.; // in case unfolding is off
1572     if (IsRecalibrationOn()) 
1573     {
1574       recalFactor = GetEMCALChannelRecalibrationFactor(iSupMod,ieta,iphi);
1575     }
1576     eCell  = cells->GetCellAmplitude(cluster->GetCellAbsId(iDigit))*fraction*recalFactor;
1577     
1578     if(cluster->E() > 0 && eCell > 0)
1579     {
1580       w  = GetCellWeight(eCell,cluster->E());
1581       
1582       etai=(Double_t)ieta;
1583       phii=(Double_t)iphi;              
1584       if(w > 0.0) 
1585       { 
1586         disp +=  w *((etai-etaMean)*(etai-etaMean)+(phii-phiMean)*(phii-phiMean)); 
1587         dEta +=  w * (etai-etaMean)*(etai-etaMean) ; 
1588         dPhi +=  w * (phii-phiMean)*(phii-phiMean) ; 
1589       }
1590     }
1591     else
1592       AliError(Form("Wrong energy %f and/or amplitude %f\n", eCell, cluster->E()));
1593   }// cell loop
1594   
1595   //Normalize to the weigth and set shower shape parameters
1596   if (wtot > 0 && nstat > 1) 
1597   {
1598     disp    /= wtot ;
1599     dEta    /= wtot ;
1600     dPhi    /= wtot ;
1601     sEta    /= wtot ;
1602     sPhi    /= wtot ;
1603     sEtaPhi /= wtot ;
1604     
1605     sEta    -= etaMean * etaMean ;
1606     sPhi    -= phiMean * phiMean ;
1607     sEtaPhi -= etaMean * phiMean ;
1608     
1609     l0 = (0.5 * (sEta + sPhi) + TMath::Sqrt( 0.25 * (sEta - sPhi) * (sEta - sPhi) + sEtaPhi * sEtaPhi ));
1610     l1 = (0.5 * (sEta + sPhi) - TMath::Sqrt( 0.25 * (sEta - sPhi) * (sEta - sPhi) + sEtaPhi * sEtaPhi ));
1611   }
1612   else
1613   {
1614     l0   = 0. ;
1615     l1   = 0. ;
1616     dEta = 0. ; dPhi = 0. ; disp    = 0. ;
1617     sEta = 0. ; sPhi = 0. ; sEtaPhi = 0. ;
1618   }     
1619   
1620 }
1621
1622 //____________________________________________________________________________________________
1623 void AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(const AliEMCALGeometry * geom, 
1624                                                                 AliVCaloCells* cells, 
1625                                                                 AliVCluster * cluster)
1626 {
1627   // Calculates new center of gravity in the local EMCAL-module coordinates 
1628   // and tranfers into global ALICE coordinates
1629   // Calculates Dispersion and main axis and puts them into the cluster
1630   
1631   Float_t l0   = 0., l1   = 0.;
1632   Float_t disp = 0., dEta = 0., dPhi    = 0.; 
1633   Float_t sEta = 0., sPhi = 0., sEtaPhi = 0.;
1634   
1635   AliEMCALRecoUtils::RecalculateClusterShowerShapeParameters(geom,cells,cluster,l0,l1,disp,
1636                                                              dEta, dPhi, sEta, sPhi, sEtaPhi);
1637   
1638   cluster->SetM02(l0);
1639   cluster->SetM20(l1);
1640   if(disp > 0. ) cluster->SetDispersion(TMath::Sqrt(disp)) ;
1641   
1642
1643
1644 //____________________________________________________________________________
1645 void AliEMCALRecoUtils::FindMatches(AliVEvent *event,
1646                                     TObjArray * clusterArr,  
1647                                     const AliEMCALGeometry *geom)
1648 {
1649   //This function should be called before the cluster loop
1650   //Before call this function, please recalculate the cluster positions
1651   //Given the input event, loop over all the tracks, select the closest cluster as matched with fCutR
1652   //Store matched cluster indexes and residuals
1653   
1654   fMatchedTrackIndex  ->Reset();
1655   fMatchedClusterIndex->Reset();
1656   fResidualPhi->Reset();
1657   fResidualEta->Reset();
1658   
1659   fMatchedTrackIndex  ->Set(1000);
1660   fMatchedClusterIndex->Set(1000);
1661   fResidualPhi->Set(1000);
1662   fResidualEta->Set(1000);
1663   
1664   AliESDEvent* esdevent = dynamic_cast<AliESDEvent*> (event);
1665   AliAODEvent* aodevent = dynamic_cast<AliAODEvent*> (event);
1666   
1667   // init the magnetic field if not already on
1668   if(!TGeoGlobalMagField::Instance()->GetField())
1669   {
1670     AliInfo("Init the magnetic field\n");
1671     if     (esdevent) 
1672     {
1673       esdevent->InitMagneticField();
1674     }
1675     else if(aodevent)
1676     {
1677       Double_t curSol = 30000*aodevent->GetMagneticField()/5.00668;
1678       Double_t curDip = 6000 *aodevent->GetMuonMagFieldScale();
1679       AliMagF *field  = AliMagF::CreateFieldMap(curSol,curDip);
1680       TGeoGlobalMagField::Instance()->SetField(field);
1681     }
1682     else
1683     {
1684       AliInfo("Mag Field not initialized, null esd/aod evetn pointers");
1685     }
1686     
1687   } // Init mag field
1688   
1689   TObjArray *clusterArray = 0x0;
1690   if(!clusterArr)
1691     {
1692       clusterArray = new TObjArray(event->GetNumberOfCaloClusters());
1693       for(Int_t icl=0; icl<event->GetNumberOfCaloClusters(); icl++)
1694         {
1695           AliVCluster *cluster = (AliVCluster*) event->GetCaloCluster(icl);
1696           if(geom && !IsGoodCluster(cluster,geom,(AliVCaloCells*)event->GetEMCALCells())) continue;
1697           clusterArray->AddAt(cluster,icl);
1698         }
1699     }
1700   
1701   Int_t    matched=0;
1702   Double_t cv[21];
1703   for (Int_t i=0; i<21;i++) cv[i]=0;
1704   for(Int_t itr=0; itr<event->GetNumberOfTracks(); itr++)
1705   {
1706     AliExternalTrackParam *trackParam = 0;
1707     
1708     //If the input event is ESD, the starting point for extrapolation is TPCOut, if available, or TPCInner 
1709     AliESDtrack *esdTrack = 0;
1710     AliAODTrack *aodTrack = 0;
1711     if(esdevent)
1712     {
1713       esdTrack = esdevent->GetTrack(itr);
1714       if(!esdTrack) continue;
1715       if(!IsAccepted(esdTrack)) continue;
1716       if(esdTrack->Pt()<fCutMinTrackPt) continue;
1717       Double_t phi = esdTrack->Phi()*TMath::RadToDeg();
1718       if(TMath::Abs(esdTrack->Eta())>0.8 || phi <= 20 || phi >= 240 ) continue;
1719       trackParam =  const_cast<AliExternalTrackParam*>(esdTrack->GetInnerParam());
1720     }
1721     
1722     //If the input event is AOD, the starting point for extrapolation is at vertex
1723     //AOD tracks are selected according to its filterbit.
1724     else if(aodevent)
1725     {
1726       aodTrack = aodevent->GetTrack(itr);
1727       if(!aodTrack) continue;
1728       if(!aodTrack->TestFilterMask(fAODFilterMask)) continue; //Select AOD tracks that fulfill GetStandardITSTPCTrackCuts2010()
1729       if(aodTrack->Pt()<fCutMinTrackPt) continue;
1730       Double_t phi = aodTrack->Phi()*TMath::RadToDeg();
1731       if(TMath::Abs(aodTrack->Eta())>0.8 || phi <= 20 || phi >= 240 ) continue;
1732       Double_t pos[3],mom[3];
1733       aodTrack->GetXYZ(pos);
1734       aodTrack->GetPxPyPz(mom);
1735       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()));
1736       trackParam= new AliExternalTrackParam(pos,mom,cv,aodTrack->Charge());
1737     }
1738     
1739     //Return if the input data is not "AOD" or "ESD"
1740     else
1741     {
1742       printf("Wrong input data type! Should be \"AOD\" or \"ESD\"\n");
1743       if(clusterArray)
1744         {
1745           clusterArray->Clear();
1746           delete clusterArray;
1747         }
1748       return;
1749     }
1750     
1751     if(!trackParam) continue;
1752
1753     //Extrapolate the track to EMCal surface
1754     AliExternalTrackParam emcalParam(*trackParam);
1755     Float_t eta, phi;
1756     if(!ExtrapolateTrackToEMCalSurface(&emcalParam, 430., fMass, fStepSurface, eta, phi)) 
1757       {
1758         if(aodevent && trackParam) delete trackParam;
1759         continue;
1760       }
1761
1762 //    if(esdevent)
1763 //      {
1764 //      esdTrack->SetOuterParam(&emcalParam,AliExternalTrackParam::kMultSec);
1765 //      }
1766
1767     if(TMath::Abs(eta)>0.75 || (phi) < 70*TMath::DegToRad() || (phi) > 190*TMath::DegToRad())
1768       {
1769         if(aodevent && trackParam) delete trackParam;
1770         continue;
1771       }
1772
1773
1774     //Find matched clusters
1775     Int_t index = -1;
1776     Float_t dEta = -999, dPhi = -999;
1777     if(!clusterArr)
1778       {
1779         index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArray, dEta, dPhi);  
1780       }
1781     else
1782       {
1783         index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArr, dEta, dPhi);  
1784       }  
1785     
1786     if(index>-1)
1787     {
1788       fMatchedTrackIndex   ->AddAt(itr,matched);
1789       fMatchedClusterIndex ->AddAt(index,matched);
1790       fResidualEta         ->AddAt(dEta,matched);
1791       fResidualPhi         ->AddAt(dPhi,matched);
1792       matched++;
1793     }
1794     if(aodevent && trackParam) delete trackParam;
1795   }//track loop
1796
1797   if(clusterArray)
1798     {
1799       clusterArray->Clear();
1800       delete clusterArray;
1801     }
1802   
1803   AliDebug(2,Form("Number of matched pairs = %d !\n",matched));
1804   
1805   fMatchedTrackIndex   ->Set(matched);
1806   fMatchedClusterIndex ->Set(matched);
1807   fResidualPhi         ->Set(matched);
1808   fResidualEta         ->Set(matched);
1809 }
1810
1811 //________________________________________________________________________________
1812 Int_t AliEMCALRecoUtils::FindMatchedClusterInEvent(const AliESDtrack *track, 
1813                                                    const AliVEvent *event, 
1814                                                    const AliEMCALGeometry *geom, 
1815                                                    Float_t &dEta, Float_t &dPhi)
1816 {
1817   //
1818   // This function returns the index of matched cluster to input track
1819   // Returns -1 if no match is found
1820   Int_t index = -1;
1821   Double_t phiV = track->Phi()*TMath::RadToDeg();
1822   if(TMath::Abs(track->Eta())>0.8 || phiV <= 20 || phiV >= 240 ) return index;
1823   AliExternalTrackParam *trackParam = const_cast<AliExternalTrackParam*>(track->GetInnerParam());
1824   if(!trackParam) return index;
1825   AliExternalTrackParam emcalParam(*trackParam);
1826   Float_t eta, phi;
1827   if(!ExtrapolateTrackToEMCalSurface(&emcalParam, 430., fMass, fStepSurface, eta, phi)) return index;
1828   if(TMath::Abs(eta)>0.75 || (phi) < 70*TMath::DegToRad() || (phi) > 190*TMath::DegToRad()) return index;
1829
1830   TObjArray *clusterArr = new TObjArray(event->GetNumberOfCaloClusters());
1831
1832   for(Int_t icl=0; icl<event->GetNumberOfCaloClusters(); icl++)
1833   {
1834     AliVCluster *cluster = (AliVCluster*) event->GetCaloCluster(icl);
1835     if(geom && !IsGoodCluster(cluster,geom,(AliVCaloCells*)event->GetEMCALCells())) continue;
1836     clusterArr->AddAt(cluster,icl);
1837   }
1838
1839   index = FindMatchedClusterInClusterArr(&emcalParam, &emcalParam, clusterArr, dEta, dPhi);  
1840   clusterArr->Clear();
1841   delete clusterArr;
1842   
1843   return index;
1844 }
1845
1846 //_______________________________________________________________________________________________
1847 Int_t  AliEMCALRecoUtils::FindMatchedClusterInClusterArr(const AliExternalTrackParam *emcalParam, 
1848                                                          AliExternalTrackParam *trkParam, 
1849                                                          const TObjArray * clusterArr, 
1850                                                          Float_t &dEta, Float_t &dPhi)
1851 {
1852   // Find matched cluster in array
1853   
1854   dEta=-999, dPhi=-999;
1855   Float_t dRMax = fCutR, dEtaMax=fCutEta, dPhiMax=fCutPhi;
1856   Int_t index = -1;
1857   Float_t tmpEta=-999, tmpPhi=-999;
1858
1859   Double_t exPos[3] = {0.,0.,0.};
1860   if(!emcalParam->GetXYZ(exPos)) return index;
1861
1862   Float_t clsPos[3] = {0.,0.,0.};
1863   for(Int_t icl=0; icl<clusterArr->GetEntriesFast(); icl++)
1864     {
1865       AliVCluster *cluster = dynamic_cast<AliVCluster*> (clusterArr->At(icl)) ;
1866       if(!cluster || !cluster->IsEMCAL()) continue;
1867       cluster->GetPosition(clsPos);
1868       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));
1869       if(dR > fClusterWindow) continue;
1870
1871       AliExternalTrackParam trkPamTmp (*trkParam);//Retrieve the starting point every time before the extrapolation
1872       if(!ExtrapolateTrackToCluster(&trkPamTmp, cluster, fMass, fStepCluster, tmpEta, tmpPhi)) continue;
1873       if(fCutEtaPhiSum)
1874         {
1875           Float_t tmpR=TMath::Sqrt(tmpEta*tmpEta + tmpPhi*tmpPhi);
1876           if(tmpR<dRMax)
1877             {
1878               dRMax=tmpR;
1879               dEtaMax=tmpEta;
1880               dPhiMax=tmpPhi;
1881               index=icl;
1882             }
1883         }
1884       else if(fCutEtaPhiSeparate)
1885         {
1886           if(TMath::Abs(tmpEta)<TMath::Abs(dEtaMax) && TMath::Abs(tmpPhi)<TMath::Abs(dPhiMax))
1887             {
1888               dEtaMax = tmpEta;
1889               dPhiMax = tmpPhi;
1890               index=icl;
1891             }
1892         }
1893       else
1894         {
1895           printf("Error: please specify your cut criteria\n");
1896           printf("To cut on sqrt(dEta^2+dPhi^2), use: SwitchOnCutEtaPhiSum()\n");
1897           printf("To cut on dEta and dPhi separately, use: SwitchOnCutEtaPhiSeparate()\n");
1898           return index;
1899         }
1900     }
1901
1902   dEta=dEtaMax;
1903   dPhi=dPhiMax;
1904
1905   return index;
1906 }
1907
1908 //------------------------------------------------------------------------------------
1909 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToEMCalSurface(AliExternalTrackParam *trkParam, 
1910                                                          const Double_t emcalR,
1911                                                          const Double_t mass, 
1912                                                          const Double_t step, 
1913                                                          Float_t &eta, 
1914                                                          Float_t &phi)
1915 {
1916   //Extrapolate track to EMCAL surface
1917   
1918   eta = -999, phi = -999;
1919   if(!trkParam) return kFALSE;
1920   if(!AliTrackerBase::PropagateTrackToBxByBz(trkParam, emcalR, mass, step, kTRUE, 0.8, -1)) return kFALSE;
1921   Double_t trkPos[3] = {0.,0.,0.};
1922   if(!trkParam->GetXYZ(trkPos)) return kFALSE;
1923   TVector3 trkPosVec(trkPos[0],trkPos[1],trkPos[2]);
1924   eta = trkPosVec.Eta();
1925   phi = trkPosVec.Phi();
1926   if(phi<0)
1927     phi += 2*TMath::Pi();
1928
1929   return kTRUE;
1930 }
1931
1932 //-----------------------------------------------------------------------------------
1933 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToPosition(AliExternalTrackParam *trkParam, 
1934                                                      const Float_t *clsPos, 
1935                                                      Double_t mass, 
1936                                                      Double_t step, 
1937                                                      Float_t &tmpEta, 
1938                                                      Float_t &tmpPhi)
1939 {
1940   //
1941   //Return the residual by extrapolating a track param to a global position
1942   //
1943   tmpEta = -999;
1944   tmpPhi = -999;
1945   if(!trkParam) return kFALSE;
1946   Double_t trkPos[3] = {0.,0.,0.};
1947   TVector3 vec(clsPos[0],clsPos[1],clsPos[2]);
1948   Double_t alpha =  ((int)(vec.Phi()*TMath::RadToDeg()/20)+0.5)*20*TMath::DegToRad();
1949   vec.RotateZ(-alpha); //Rotate the cluster to the local extrapolation coordinate system
1950   if(!AliTrackerBase::PropagateTrackToBxByBz(trkParam, vec.X(), mass, step,kTRUE, 0.8, -1)) return kFALSE;
1951   if(!trkParam->GetXYZ(trkPos)) return kFALSE; //Get the extrapolated global position
1952
1953   TVector3 clsPosVec(clsPos[0],clsPos[1],clsPos[2]);
1954   TVector3 trkPosVec(trkPos[0],trkPos[1],trkPos[2]);
1955
1956   // track cluster matching
1957   tmpPhi = clsPosVec.DeltaPhi(trkPosVec);    // tmpPhi is between -pi and pi
1958   tmpEta = clsPosVec.Eta()-trkPosVec.Eta();
1959
1960   return kTRUE;
1961 }
1962
1963 //----------------------------------------------------------------------------------
1964 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToCluster(AliExternalTrackParam *trkParam, 
1965                                                     const AliVCluster *cluster, 
1966                                                     const Double_t mass, 
1967                                                     const Double_t step, 
1968                                                     Float_t &tmpEta, 
1969                                                     Float_t &tmpPhi)
1970 {
1971   //
1972   //Return the residual by extrapolating a track param to a cluster
1973   //
1974   tmpEta = -999;
1975   tmpPhi = -999;
1976   if(!cluster || !trkParam) return kFALSE;
1977
1978   Float_t clsPos[3] = {0.,0.,0.};
1979   cluster->GetPosition(clsPos);
1980
1981   return ExtrapolateTrackToPosition(trkParam, clsPos, mass, step, tmpEta, tmpPhi);
1982 }
1983
1984 //---------------------------------------------------------------------------------
1985 Bool_t AliEMCALRecoUtils::ExtrapolateTrackToCluster(AliExternalTrackParam *trkParam, 
1986                                                     const AliVCluster *cluster, 
1987                                                     Float_t &tmpEta, 
1988                                                     Float_t &tmpPhi)
1989 {
1990   //
1991   //Return the residual by extrapolating a track param to a clusterfStepCluster
1992   //
1993
1994   return ExtrapolateTrackToCluster(trkParam, cluster, fMass, fStepCluster, tmpEta, tmpPhi);
1995 }
1996
1997 //_______________________________________________________________________
1998 void AliEMCALRecoUtils::GetMatchedResiduals(const Int_t clsIndex, 
1999                                             Float_t &dEta, Float_t &dPhi)
2000 {
2001   //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2002   //Get the residuals dEta and dPhi for this cluster to the closest track
2003   //Works with ESDs and AODs
2004
2005   if( FindMatchedPosForCluster(clsIndex) >= 999 )
2006   {
2007     AliDebug(2,"No matched tracks found!\n");
2008     dEta=999.;
2009     dPhi=999.;
2010     return;
2011   }
2012   dEta = fResidualEta->At(FindMatchedPosForCluster(clsIndex));
2013   dPhi = fResidualPhi->At(FindMatchedPosForCluster(clsIndex));
2014 }
2015
2016 //______________________________________________________________________________________________
2017 void AliEMCALRecoUtils::GetMatchedClusterResiduals(Int_t trkIndex, Float_t &dEta, Float_t &dPhi)
2018 {
2019   //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2020   //Get the residuals dEta and dPhi for this track to the closest cluster
2021   //Works with ESDs and AODs
2022
2023   if( FindMatchedPosForTrack(trkIndex) >= 999 )
2024   {
2025     AliDebug(2,"No matched cluster found!\n");
2026     dEta=999.;
2027     dPhi=999.;
2028     return;
2029   }
2030   dEta = fResidualEta->At(FindMatchedPosForTrack(trkIndex));
2031   dPhi = fResidualPhi->At(FindMatchedPosForTrack(trkIndex));
2032 }
2033
2034 //__________________________________________________________
2035 Int_t AliEMCALRecoUtils::GetMatchedTrackIndex(Int_t clsIndex)
2036 {
2037   //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2038   //Get the index of matched track to this cluster
2039   //Works with ESDs and AODs
2040   
2041   if(IsClusterMatched(clsIndex))
2042     return fMatchedTrackIndex->At(FindMatchedPosForCluster(clsIndex));
2043   else 
2044     return -1; 
2045 }
2046
2047 //__________________________________________________________
2048 Int_t AliEMCALRecoUtils::GetMatchedClusterIndex(Int_t trkIndex)
2049 {
2050   //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2051   //Get the index of matched cluster to this track
2052   //Works with ESDs and AODs
2053   
2054   if(IsTrackMatched(trkIndex))
2055     return fMatchedClusterIndex->At(FindMatchedPosForTrack(trkIndex));
2056   else 
2057     return -1; 
2058 }
2059
2060 //______________________________________________________________
2061 Bool_t AliEMCALRecoUtils::IsClusterMatched(Int_t clsIndex) const
2062 {
2063   //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2064   //Returns if the cluster has a match
2065   if(FindMatchedPosForCluster(clsIndex) < 999) 
2066     return kTRUE;
2067   else
2068     return kFALSE;
2069 }
2070
2071 //____________________________________________________________
2072 Bool_t AliEMCALRecoUtils::IsTrackMatched(Int_t trkIndex) const 
2073 {
2074   //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2075   //Returns if the track has a match
2076   if(FindMatchedPosForTrack(trkIndex) < 999) 
2077     return kTRUE;
2078   else
2079     return kFALSE;
2080 }
2081
2082 //______________________________________________________________________
2083 UInt_t AliEMCALRecoUtils::FindMatchedPosForCluster(Int_t clsIndex) const
2084 {
2085   //Given a cluster index as in AliESDEvent::GetCaloCluster(clsIndex)
2086   //Returns the position of the match in the fMatchedClusterIndex array
2087   Float_t tmpR = fCutR;
2088   UInt_t pos = 999;
2089   
2090   for(Int_t i=0; i<fMatchedClusterIndex->GetSize(); i++) 
2091   {
2092     if(fMatchedClusterIndex->At(i)==clsIndex) 
2093     {
2094       Float_t r = TMath::Sqrt(fResidualEta->At(i)*fResidualEta->At(i) + fResidualPhi->At(i)*fResidualPhi->At(i));
2095       if(r<tmpR) 
2096       {
2097         pos=i;
2098         tmpR=r;
2099         AliDebug(3,Form("Matched cluster index: index: %d, dEta: %2.4f, dPhi: %2.4f.\n",
2100                         fMatchedClusterIndex->At(i),fResidualEta->At(i),fResidualPhi->At(i)));
2101       }
2102     }
2103   }
2104   return pos;
2105 }
2106
2107 //____________________________________________________________________
2108 UInt_t AliEMCALRecoUtils::FindMatchedPosForTrack(Int_t trkIndex) const
2109 {
2110   //Given a track index as in AliESDEvent::GetTrack(trkIndex)
2111   //Returns the position of the match in the fMatchedTrackIndex array
2112   Float_t tmpR = fCutR;
2113   UInt_t pos = 999;
2114   
2115   for(Int_t i=0; i<fMatchedTrackIndex->GetSize(); i++) 
2116   {
2117     if(fMatchedTrackIndex->At(i)==trkIndex) 
2118     {
2119       Float_t r = TMath::Sqrt(fResidualEta->At(i)*fResidualEta->At(i) + fResidualPhi->At(i)*fResidualPhi->At(i));
2120       if(r<tmpR) 
2121       {
2122         pos=i;
2123         tmpR=r;
2124         AliDebug(3,Form("Matched track index: index: %d, dEta: %2.4f, dPhi: %2.4f.\n",
2125                         fMatchedTrackIndex->At(i),fResidualEta->At(i),fResidualPhi->At(i)));
2126       }
2127     }
2128   }
2129   return pos;
2130 }
2131
2132 //__________________________________________________________________________
2133 Bool_t AliEMCALRecoUtils::IsGoodCluster(AliVCluster *cluster, 
2134                                         const AliEMCALGeometry *geom, 
2135                                         AliVCaloCells* cells,const Int_t bc)
2136 {
2137   // check if the cluster survives some quality cut
2138   //
2139   //
2140   Bool_t isGood=kTRUE;
2141
2142   if(!cluster || !cluster->IsEMCAL())              return kFALSE;
2143   
2144   if(ClusterContainsBadChannel(geom,cluster->GetCellsAbsId(),cluster->GetNCells())) return kFALSE;
2145   
2146   if(!CheckCellFiducialRegion(geom,cluster,cells)) return kFALSE;
2147   
2148   if(IsExoticCluster(cluster, cells,bc))           return kFALSE;
2149
2150   return isGood;
2151 }
2152
2153 //__________________________________________________________
2154 Bool_t AliEMCALRecoUtils::IsAccepted(AliESDtrack *esdTrack)
2155 {
2156   // Given a esd track, return whether the track survive all the cuts
2157
2158   // The different quality parameter are first
2159   // retrieved from the track. then it is found out what cuts the
2160   // track did not survive and finally the cuts are imposed.
2161
2162   UInt_t status = esdTrack->GetStatus();
2163
2164   Int_t nClustersITS = esdTrack->GetITSclusters(0);
2165   Int_t nClustersTPC = esdTrack->GetTPCclusters(0);
2166
2167   Float_t chi2PerClusterITS = -1;
2168   Float_t chi2PerClusterTPC = -1;
2169   if (nClustersITS!=0)
2170     chi2PerClusterITS = esdTrack->GetITSchi2()/Float_t(nClustersITS);
2171   if (nClustersTPC!=0) 
2172     chi2PerClusterTPC = esdTrack->GetTPCchi2()/Float_t(nClustersTPC);
2173
2174
2175   //DCA cuts
2176   if(fTrackCutsType==kGlobalCut)
2177     {
2178       Float_t maxDCAToVertexXYPtDep = 0.0182 + 0.0350/TMath::Power(esdTrack->Pt(),1.01); //This expression comes from AliESDtrackCuts::GetStandardITSTPCTrackCuts2010()
2179       //AliDebug(3,Form("Track pT = %f, DCAtoVertexXY = %f",esdTrack->Pt(),MaxDCAToVertexXYPtDep));
2180       SetMaxDCAToVertexXY(maxDCAToVertexXYPtDep); //Set pT dependent DCA cut to vertex in x-y plane
2181     }
2182
2183
2184   Float_t b[2];
2185   Float_t bCov[3];
2186   esdTrack->GetImpactParameters(b,bCov);
2187   if (bCov[0]<=0 || bCov[2]<=0) 
2188   {
2189     AliDebug(1, "Estimated b resolution lower or equal zero!");
2190     bCov[0]=0; bCov[2]=0;
2191   }
2192
2193   Float_t dcaToVertexXY = b[0];
2194   Float_t dcaToVertexZ = b[1];
2195   Float_t dcaToVertex = -1;
2196
2197   if (fCutDCAToVertex2D)
2198     dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY/fCutMaxDCAToVertexXY/fCutMaxDCAToVertexXY + dcaToVertexZ*dcaToVertexZ/fCutMaxDCAToVertexZ/fCutMaxDCAToVertexZ);
2199   else
2200     dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY + dcaToVertexZ*dcaToVertexZ);
2201     
2202   // cut the track?
2203   
2204   Bool_t cuts[kNCuts];
2205   for (Int_t i=0; i<kNCuts; i++) cuts[i]=kFALSE;
2206   
2207   // track quality cuts
2208   if (fCutRequireTPCRefit && (status&AliESDtrack::kTPCrefit)==0)
2209     cuts[0]=kTRUE;
2210   if (fCutRequireITSRefit && (status&AliESDtrack::kITSrefit)==0)
2211     cuts[1]=kTRUE;
2212   if (nClustersTPC<fCutMinNClusterTPC)
2213     cuts[2]=kTRUE;
2214   if (nClustersITS<fCutMinNClusterITS) 
2215     cuts[3]=kTRUE;
2216   if (chi2PerClusterTPC>fCutMaxChi2PerClusterTPC) 
2217     cuts[4]=kTRUE; 
2218   if (chi2PerClusterITS>fCutMaxChi2PerClusterITS) 
2219     cuts[5]=kTRUE;  
2220   if (!fCutAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0)
2221     cuts[6]=kTRUE;
2222   if (fCutDCAToVertex2D && dcaToVertex > 1)
2223     cuts[7] = kTRUE;
2224   if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexXY) > fCutMaxDCAToVertexXY)
2225     cuts[8] = kTRUE;
2226   if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexZ) > fCutMaxDCAToVertexZ)
2227     cuts[9] = kTRUE;
2228
2229   if(fTrackCutsType==kGlobalCut)
2230     {
2231       //Require at least one SPD point + anything else in ITS
2232       if( (esdTrack->HasPointOnITSLayer(0) || esdTrack->HasPointOnITSLayer(1)) == kFALSE)
2233         cuts[10] = kTRUE;
2234     }
2235
2236   Bool_t cut=kFALSE;
2237   for (Int_t i=0; i<kNCuts; i++)
2238     if (cuts[i]) { cut = kTRUE ; }
2239
2240     // cut the track
2241   if (cut) 
2242     return kFALSE;
2243   else 
2244     return kTRUE;
2245 }
2246
2247 //_____________________________________
2248 void AliEMCALRecoUtils::InitTrackCuts()
2249 {
2250   //Intilize the track cut criteria
2251   //By default these cuts are set according to AliESDtrackCuts::GetStandardTPCOnlyTrackCuts()
2252   //Also you can customize the cuts using the setters
2253   
2254   switch (fTrackCutsType)
2255   {
2256     case kTPCOnlyCut:
2257     {
2258       AliInfo(Form("Track cuts for matching: GetStandardTPCOnlyTrackCuts()"));
2259       //TPC
2260       SetMinNClustersTPC(70);
2261       SetMaxChi2PerClusterTPC(4);
2262       SetAcceptKinkDaughters(kFALSE);
2263       SetRequireTPCRefit(kFALSE);
2264       
2265       //ITS
2266       SetRequireITSRefit(kFALSE);
2267       SetMaxDCAToVertexZ(3.2);
2268       SetMaxDCAToVertexXY(2.4);
2269       SetDCAToVertex2D(kTRUE);
2270       
2271       break;
2272     }
2273       
2274     case kGlobalCut:
2275     {
2276       AliInfo(Form("Track cuts for matching: GetStandardITSTPCTrackCuts2010(kTURE)"));
2277       //TPC
2278       SetMinNClustersTPC(70);
2279       SetMaxChi2PerClusterTPC(4);
2280       SetAcceptKinkDaughters(kFALSE);
2281       SetRequireTPCRefit(kTRUE);
2282       
2283       //ITS
2284       SetRequireITSRefit(kTRUE);
2285       SetMaxDCAToVertexZ(2);
2286       SetMaxDCAToVertexXY();
2287       SetDCAToVertex2D(kFALSE);
2288       
2289       break;
2290     }
2291       
2292     case kLooseCut:
2293     {
2294       AliInfo(Form("Track cuts for matching: Loose cut w/o DCA cut"));
2295       SetMinNClustersTPC(50);
2296       SetAcceptKinkDaughters(kTRUE);
2297       
2298       break;
2299     }
2300   }
2301 }
2302
2303
2304 //________________________________________________________________________
2305 void AliEMCALRecoUtils::SetClusterMatchedToTrack(const AliVEvent *event)
2306 {
2307   // Checks if tracks are matched to EMC clusters and set the matched EMCAL cluster index to ESD track. 
2308
2309   Int_t nTracks = event->GetNumberOfTracks();
2310   for (Int_t iTrack = 0; iTrack < nTracks; ++iTrack) 
2311   {
2312     AliVTrack* track = dynamic_cast<AliVTrack*>(event->GetTrack(iTrack));
2313     if (!track) 
2314     {
2315       AliWarning(Form("Could not receive track %d", iTrack));
2316       continue;
2317     }
2318     
2319     Int_t matchClusIndex = GetMatchedClusterIndex(iTrack);                 
2320     track->SetEMCALcluster(matchClusIndex); //sets -1 if track not matched within residual
2321     /*the following can be done better if AliVTrack::SetStatus will be there. Patch pending with Andreas/Peter*/
2322     AliESDtrack* esdtrack = dynamic_cast<AliESDtrack*>(track);
2323     if (esdtrack) { 
2324       if(matchClusIndex != -1) 
2325         esdtrack->SetStatus(AliESDtrack::kEMCALmatch);
2326       else
2327         esdtrack->ResetStatus(AliESDtrack::kEMCALmatch);
2328     } else {
2329       AliAODTrack* aodtrack = dynamic_cast<AliAODTrack*>(track);
2330       if(matchClusIndex != -1) 
2331         aodtrack->SetStatus(AliESDtrack::kEMCALmatch);
2332       else
2333         aodtrack->ResetStatus(AliESDtrack::kEMCALmatch);
2334     }
2335
2336   }
2337   AliDebug(2,"Track matched to closest cluster");       
2338 }
2339
2340 //_________________________________________________________________________
2341 void AliEMCALRecoUtils::SetTracksMatchedToCluster(const AliVEvent *event)
2342 {
2343   // Checks if EMC clusters are matched to ESD track.
2344   // Adds track indexes of all the tracks matched to a cluster withing residuals in ESDCalocluster.
2345   
2346   for (Int_t iClus=0; iClus < event->GetNumberOfCaloClusters(); ++iClus) 
2347   {
2348     AliVCluster *cluster = event->GetCaloCluster(iClus);
2349     if (!cluster->IsEMCAL()) 
2350       continue;
2351     
2352     Int_t nTracks = event->GetNumberOfTracks();
2353     TArrayI arrayTrackMatched(nTracks);
2354     
2355     // Get the closest track matched to the cluster
2356     Int_t nMatched = 0;
2357     Int_t matchTrackIndex = GetMatchedTrackIndex(iClus);
2358     if (matchTrackIndex != -1) 
2359     {
2360       arrayTrackMatched[nMatched] = matchTrackIndex;
2361       nMatched++;
2362     }
2363     
2364     // Get all other tracks matched to the cluster
2365     for(Int_t iTrk=0; iTrk<nTracks; ++iTrk) 
2366     {
2367       AliVTrack* track = dynamic_cast<AliVTrack*>(event->GetTrack(iTrk));
2368       if(iTrk == matchTrackIndex) continue;
2369       if(track->GetEMCALcluster() == iClus)
2370       {
2371         arrayTrackMatched[nMatched] = iTrk;
2372         ++nMatched;
2373       }
2374     }
2375     
2376     //printf("Tender::SetTracksMatchedToCluster - cluster E %f, N matches %d, first match %d\n",cluster->E(),nMatched,arrayTrackMatched[0]);
2377     
2378     arrayTrackMatched.Set(nMatched);
2379     AliESDCaloCluster *esdcluster = dynamic_cast<AliESDCaloCluster*>(cluster);
2380     if (esdcluster) 
2381       esdcluster->AddTracksMatched(arrayTrackMatched);
2382     else if (nMatched>0) {
2383       AliAODCaloCluster *aodcluster = dynamic_cast<AliAODCaloCluster*>(cluster);
2384       if (aodcluster)
2385         aodcluster->AddTrackMatched(event->GetTrack(arrayTrackMatched.At(0)));
2386     }
2387     
2388     Float_t eta= -999, phi = -999;
2389     if (matchTrackIndex != -1) 
2390       GetMatchedResiduals(iClus, eta, phi);
2391     cluster->SetTrackDistance(phi, eta);
2392   }
2393   
2394   AliDebug(2,"Cluster matched to tracks");      
2395 }
2396
2397
2398 //___________________________________________________
2399 void AliEMCALRecoUtils::Print(const Option_t *) const 
2400 {
2401   // Print Parameters
2402   
2403   printf("AliEMCALRecoUtils Settings: \n");
2404   printf("Misalignment shifts\n");
2405   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, 
2406                                   fMisalTransShift[i*3],fMisalTransShift[i*3+1],fMisalTransShift[i*3+2],
2407                                   fMisalRotShift[i*3],  fMisalRotShift[i*3+1],  fMisalRotShift[i*3+2]   );
2408   printf("Non linearity function %d, parameters:\n", fNonLinearityFunction);
2409   for(Int_t i=0; i<6; i++) printf("param[%d]=%f\n",i, fNonLinearityParams[i]);
2410   
2411   printf("Position Recalculation option %d, Particle Type %d, fW0 %2.2f, Recalibrate Data %d \n",fPosAlgo,fParticleType,fW0, fRecalibration);
2412
2413   printf("Matching criteria: ");
2414   if(fCutEtaPhiSum)
2415     {
2416       printf("sqrt(dEta^2+dPhi^2)<%4.3f\n",fCutR);
2417     }
2418   else if(fCutEtaPhiSeparate)
2419     {
2420       printf("dEta<%4.3f, dPhi<%4.3f\n",fCutEta,fCutPhi);
2421     }
2422   else
2423     {
2424       printf("Error\n");
2425       printf("please specify your cut criteria\n");
2426       printf("To cut on sqrt(dEta^2+dPhi^2), use: SwitchOnCutEtaPhiSum()\n");
2427       printf("To cut on dEta and dPhi separately, use: SwitchOnCutEtaPhiSeparate()\n");
2428     }
2429
2430   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);
2431   printf("Cluster selection window: dR < %2.0f\n",fClusterWindow);
2432
2433   printf("Track cuts: \n");
2434   printf("Minimum track pT: %1.2f\n",fCutMinTrackPt);
2435   printf("AOD track selection mask: %d\n",fAODFilterMask);
2436   printf("TPCRefit = %d, ITSRefit = %d\n",fCutRequireTPCRefit,fCutRequireITSRefit);
2437   printf("AcceptKinks = %d\n",fCutAcceptKinkDaughters);
2438   printf("MinNCulsterTPC = %d, MinNClusterITS = %d\n",fCutMinNClusterTPC,fCutMinNClusterITS);
2439   printf("MaxChi2TPC = %2.2f, MaxChi2ITS = %2.2f\n",fCutMaxChi2PerClusterTPC,fCutMaxChi2PerClusterITS);
2440   printf("DCSToVertex2D = %d, MaxDCAToVertexXY = %2.2f, MaxDCAToVertexZ = %2.2f\n",fCutDCAToVertex2D,fCutMaxDCAToVertexXY,fCutMaxDCAToVertexZ);
2441 }
2442
2443 //_________________________________________________________________
2444 void AliEMCALRecoUtils::SetRunDependentCorrections(Int_t runnumber)
2445 {
2446   //Get EMCAL time dependent corrections from file and put them in the recalibration histograms
2447   //Do it only once and only if it is requested
2448   
2449   if(!fUseRunCorrectionFactors) return;
2450   if(fRunCorrectionFactorsSet)  return;
2451   
2452   AliInfo(Form("AliEMCALRecoUtils::GetRunDependentCorrections() - Get Correction Factors for Run number %d\n",runnumber));
2453  
2454   AliEMCALCalibTimeDepCorrection  *corr =  new AliEMCALCalibTimeDepCorrection();
2455   corr->ReadRootInfo(Form("CorrectionFiles/Run%d_Correction.root",runnumber));
2456   
2457   SwitchOnRecalibration();
2458   
2459   AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
2460   
2461   for(Int_t ism = 0; ism < geom->GetNumberOfSuperModules(); ism++)
2462   {
2463     for(Int_t icol = 0; icol < 48; icol++)
2464     {
2465       for(Int_t irow = 0; irow < 24; irow++)
2466       {
2467         Float_t orgRecalFactor = GetEMCALChannelRecalibrationFactors(ism)->GetBinContent(icol,irow);
2468         Float_t newRecalFactor = orgRecalFactor*corr->GetCorrection(ism, icol,irow,0);
2469         GetEMCALChannelRecalibrationFactors(ism)->SetBinContent(icol,irow,newRecalFactor);
2470         //printf("ism %d, icol %d, irow %d, corrections : org %f, time dep %f, final %f (org*time %f)\n",ism, icol, irow, 
2471         //        orgRecalFactor, corr->GetCorrection(ism, icol,irow,0),
2472         //       (GetEMCALChannelRecalibrationFactors(ism))->GetBinContent(icol,irow),newRecalFactor);
2473       }
2474     }
2475   }
2476   
2477   fRunCorrectionFactorsSet = kTRUE;
2478   
2479 }