ROOT logo
#ifndef ALICALORIMETERUTILS_H
#define ALICALORIMETERUTILS_H
/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
 * See cxx source for full Copyright notice     */

//_________________________________________________________________________
// Class utility for Calorimeter specific selection methods                ///
//
//
//
//-- Author: Gustavo Conesa (LPSC-Grenoble) 
//////////////////////////////////////////////////////////////////////////////

// --- ROOT system ---
#include <TObject.h> 
#include <TString.h>
#include <TObjArray.h>
class TArrayF;  
#include <TH2I.h>
#include <TGeoMatrix.h>

//--- ANALYSIS system ---
class AliVEvent;
class AliVTrack;
class AliAODPWG4Particle;
class AliAODCaloCluster;
class AliVCaloCells;
class AliPHOSGeoUtils;
class AliEMCALGeometry;
class AliAODMCParticle;
class TParticle;

#include "AliEMCALRecoUtils.h"

class AliCalorimeterUtils : public TObject {

 public:   
  AliCalorimeterUtils() ; // ctor
  virtual ~AliCalorimeterUtils() ;//virtual dtor
  
  virtual void  InitParameters();
  virtual void  Print(const Option_t * opt)          const ;

  virtual Int_t GetDebug()                           const { return fDebug                ; }
  virtual void  SetDebug(Int_t d)                          { fDebug = d                   ; }
	
  //virtual void Init();
	
  // Cluster contents
  
  Bool_t        AreNeighbours(Int_t calo, Int_t absId1, Int_t absId2) const ;

  Bool_t        IsClusterSharedByTwoSuperModules(const AliEMCALGeometry * geom,
                                                 AliVCluster* cluster);
  
  Int_t         GetNumberOfLocalMaxima(AliVCluster* cluster, AliVCaloCells* cells)  ;
  
  Int_t         GetNumberOfLocalMaxima(AliVCluster* cluster, AliVCaloCells* cells,
                                       Int_t *absIdList,     Float_t *maxEList)  ;
  
  Float_t       GetLocalMaximaCutE()                 const { return fLocMaxCutE           ; }
  void          SetLocalMaximaCutE(Float_t cut)            { fLocMaxCutE     = cut        ; }
  
  Float_t       GetLocalMaximaCutEDiff()             const { return fLocMaxCutEDiff       ; }
  void          SetLocalMaximaCutEDiff(Float_t c)          { fLocMaxCutEDiff = c          ; }
  
  Int_t         GetMaxEnergyCell(AliVCaloCells* cells, const AliVCluster* clu, Float_t & fraction) const ;
  
  void          SplitEnergy(Int_t absId1, Int_t absId2, AliVCluster *cluster, AliVCaloCells* cells,
                            //Float_t & e1, Float_t & e2,
                            AliAODCaloCluster *cluster1, AliAODCaloCluster *cluster2,
                            Int_t nMax, Int_t eventNumber = 0);//, Int_t *absIdList, Float_t *maxEList,
  
  void          SwitchOnClusterPlot()                      { fPlotCluster = kTRUE         ; }
  void          SwitchOffClusterPlot()                     { fPlotCluster = kFALSE        ; }

  Float_t       GetMCECellClusFracCorrection(Float_t eCell, Float_t eCluster) const ;
  void          SetMCECellClusFracCorrectionParamters(Int_t i, Float_t param) { if(i<4) fMCECellClusFracCorrParam[i] = param; }
  
  Bool_t        IsMCECellClusFracCorrectionOn()   const    { return fMCECellClusFracCorrOn   ; }
  void          SwitchOnMCECellClusFracCorrection()        { fMCECellClusFracCorrOn = kTRUE  ; }
  void          SwitchOffMCECellClusFracCorrection()       { fMCECellClusFracCorrOn = kFALSE ; }

  
  //Calorimeters Geometry Methods
  AliEMCALGeometry * GetEMCALGeometry()              const { return fEMCALGeo             ; }
  TString       EMCALGeometryName()                  const { return fEMCALGeoName         ; }  
  void          SetEMCALGeometryName(TString name)         { fEMCALGeoName = name         ; }
  void          InitEMCALGeometry(Int_t runnumber = 180000) ; 
  Bool_t        IsEMCALGeoMatrixSet()                const { return fEMCALGeoMatrixSet    ; }
	
  AliPHOSGeoUtils * GetPHOSGeometry()                const { return fPHOSGeo              ; }	
  TString       PHOSGeometryName()                   const { return fPHOSGeoName          ; }  
  void          SetPHOSGeometryName(TString name)          { fPHOSGeoName = name          ; }
  void          InitPHOSGeometry(Int_t runnumber = 180000) ; 
  Bool_t        IsPHOSGeoMatrixSet()                 const { return fPHOSGeoMatrixSet     ; }

  void          AccessGeometry(AliVEvent* inputEvent) ;
	
  void          SetImportGeometryFromFile(Bool_t import,
                                          TString path = ""){
                                                             fImportGeometryFromFile = import    ;
                                                             fImportGeometryFilePath = path      ; } // EMCAL
  
  Bool_t        IsMCParticleInCalorimeterAcceptance(Int_t calo, TParticle* particle);
  Bool_t        IsMCParticleInCalorimeterAcceptance(Int_t calo, AliAODMCParticle* particle);
  Bool_t        IsMCParticleInCalorimeterAcceptance(Int_t calo, Float_t eta, Float_t theta, Float_t phi, Int_t & absID);
  
  void          SwitchOnLoadOwnEMCALGeometryMatrices()     { fLoadEMCALMatrices = kTRUE   ; }
  void          SwitchOffLoadOwnEMCALGeometryMatrices()    { fLoadEMCALMatrices = kFALSE  ; }
  void          SetEMCALGeometryMatrixInSM(TGeoHMatrix* m, Int_t i) { fEMCALMatrix[i] = m ; }
  
  void          SwitchOnLoadOwnPHOSGeometryMatrices()      { fLoadPHOSMatrices = kTRUE    ; }
  void          SwitchOffLoadOwnPHOSGeometryMatrices()     { fLoadPHOSMatrices = kFALSE   ; }
  void          SetPHOSGeometryMatrixInSM(TGeoHMatrix* m, Int_t i) { fPHOSMatrix[i] = m   ; }
  
  // Bad channels
  Bool_t        IsBadChannelsRemovalSwitchedOn()     const { return fRemoveBadChannels                                ; }
  void          SwitchOnBadChannelsRemoval ()              { fRemoveBadChannels = kTRUE   ; 
                                                             fEMCALRecoUtils->SwitchOnBadChannelsRemoval(); 
                                                             if(!fPHOSBadChannelMap) InitPHOSBadChannelStatusMap()    ; }
  void          SwitchOffBadChannelsRemoval()              { fRemoveBadChannels = kFALSE ; 
                                                             fEMCALRecoUtils->SwitchOffBadChannelsRemoval()           ; }
  
  Bool_t        IsDistanceToBadChannelRecalculated() const { return  IsDistanceToBadChannelRecalculated()             ; }
  void          SwitchOnDistToBadChannelRecalculation ()   { fEMCALRecoUtils->SwitchOnDistToBadChannelRecalculation() ; }
  void          SwitchOffDistToBadChannelRecalculation()   { fEMCALRecoUtils->SwitchOffDistToBadChannelRecalculation(); }
  
  void          InitPHOSBadChannelStatusMap () ;

  Int_t         GetEMCALChannelStatus(Int_t iSM , Int_t iCol, Int_t iRow) const { 
                  return fEMCALRecoUtils->GetEMCALChannelStatus(iSM,iCol,iRow); }//Channel is ok by default

  Int_t         GetPHOSChannelStatus (Int_t imod, Int_t iCol, Int_t iRow) const { 
                  if(fPHOSBadChannelMap) return (Int_t) ((TH2I*)fPHOSBadChannelMap->At(imod))->GetBinContent(iCol,iRow); 
                  else return 0 ; }//Channel is ok by default
  
  void          SetEMCALChannelStatus(Int_t iSM , Int_t iCol, Int_t iRow, Double_t c = 1) { 
                  fEMCALRecoUtils->SetEMCALChannelStatus(iSM,iCol,iRow,c) ; }
  
  void          SetPHOSChannelStatus (Int_t imod, Int_t iCol, Int_t iRow, Double_t c = 1) {
                  if(!fPHOSBadChannelMap) InitPHOSBadChannelStatusMap() ; 
                  ((TH2I*)fPHOSBadChannelMap->At(imod))->SetBinContent(iCol,iRow,c) ; }
    
  void          SetEMCALChannelStatusMap(Int_t iSM , TH2I* h) { fEMCALRecoUtils->SetEMCALChannelStatusMap(iSM,h)      ; }
  void          SetPHOSChannelStatusMap(Int_t imod , TH2I* h) { fPHOSBadChannelMap ->AddAt(h,imod)                    ; }
  
  TH2I *        GetEMCALChannelStatusMap(Int_t iSM)  const { return fEMCALRecoUtils->GetEMCALChannelStatusMap(iSM)    ; }
  TH2I *        GetPHOSChannelStatusMap(Int_t imod)  const { return (TH2I*)fPHOSBadChannelMap->At(imod)               ; }

  void          SetEMCALChannelStatusMap(TObjArray *map)   { fEMCALRecoUtils->SetEMCALChannelStatusMap(map)           ; }
  void          SetPHOSChannelStatusMap (TObjArray *map)   { fPHOSBadChannelMap  = map                                ; }
	
  Bool_t        ClusterContainsBadChannel(Int_t calo,UShort_t* cellList, Int_t nCells);
  Bool_t        ClusterContainsBadChannel(TString /*calo*/,UShort_t* /*cellList*/, Int_t /*nCells*/); // Stupid thing to do but just to avoid compilation break in AliTrackComparisonESD while I find the authors
	
  // Mask clusters in front of frame, EMCAL only
  Int_t         GetNMaskCellColumns()                const { return fNMaskCellColumns;}
  void          SetNMaskCellColumns(Int_t n) {
                  if(n > fNMaskCellColumns) { delete [] fMaskCellColumns ; fMaskCellColumns = new Int_t[n] ; }
                  fNMaskCellColumns = n                                                                               ; }
  void          SetMaskCellColumn(Int_t ipos, Int_t icol) { 
                  if(ipos < fNMaskCellColumns) fMaskCellColumns[ipos] = icol;
                  else printf("Not set, position larger than allocated set size first")                               ; }
  Bool_t        MaskFrameCluster(Int_t iSM, Int_t ieta) const ;
  
  
  //Calorimeter indexes information
  Int_t         GetModuleNumber(AliAODPWG4Particle * particle, AliVEvent* inputEvent) const;
  Int_t         GetModuleNumber(AliVCluster * cluster) const;
  Int_t         GetModuleNumberCellIndexes(Int_t absId, Int_t calo, Int_t & icol, Int_t & irow, Int_t &iRCU) const ;
	
  //Modules fiducial region
  Bool_t        CheckCellFiducialRegion(AliVCluster* cluster, AliVCaloCells* cells) const ;
  Bool_t        CheckCellFiducialRegion(AliVCluster* cluster, AliVCaloCells* cells, AliVEvent* /**/, Int_t /**/)
                { return CheckCellFiducialRegion(cluster, cells) ; } // Stupid thing to do but just to avoid compilation break in AliTrackComparisonESD while I find the authors
  void          SetNumberOfCellsFromPHOSBorder(Int_t n)    { fNCellsFromPHOSBorder = n                                ; }
  Int_t         GetNumberOfCellsFromPHOSBorder()     const { return fNCellsFromPHOSBorder                             ; }
  void          SetNumberOfCellsFromEMCALBorder(Int_t n)   { fEMCALRecoUtils->SetNumberOfCellsFromEMCALBorder(n)      ; }
  Int_t         GetNumberOfCellsFromEMCALBorder()    const { return fEMCALRecoUtils->GetNumberOfCellsFromEMCALBorder(); }
  void          SwitchOnNoFiducialBorderInEMCALEta0()      { fEMCALRecoUtils->SwitchOnNoFiducialBorderInEMCALEta0()   ; }
  void          SwitchOffNoFiducialBorderInEMCALEta0()     { fEMCALRecoUtils->SwitchOffNoFiducialBorderInEMCALEta0()  ; }
  Bool_t        IsEMCALNoBorderAtEta0()              const { return fEMCALRecoUtils->IsEMCALNoBorderAtEta0()          ; }
  
  // Recalibration
  Bool_t        IsRecalibrationOn()                  const { return fRecalibration                                    ; }
  void          SwitchOnRecalibration()                    { fRecalibration = kTRUE ; 
                  InitPHOSRecalibrationFactors(); fEMCALRecoUtils->SwitchOnRecalibration()                            ; }
  void          SwitchOffRecalibration()                   { fRecalibration = kFALSE;
                  fEMCALRecoUtils->SwitchOffRecalibration()                                                           ; }
	
  void          InitPHOSRecalibrationFactors () ;
	
  Float_t       GetEMCALChannelRecalibrationFactor(Int_t iSM , Int_t iCol, Int_t iRow) const { 
                  return fEMCALRecoUtils->GetEMCALChannelRecalibrationFactor(iSM , iCol, iRow)                        ; }
  
  Float_t       GetPHOSChannelRecalibrationFactor (Int_t imod, Int_t iCol, Int_t iRow) const { 
                  if(fPHOSRecalibrationFactors)
                    return (Float_t) ((TH2F*)fPHOSRecalibrationFactors->At(imod))->GetBinContent(iCol,iRow); 
                  else return 1                                                                                       ; }
  
  void          SetEMCALChannelRecalibrationFactor(Int_t iSM , Int_t iCol, Int_t iRow, Double_t c = 1) { 
                  fEMCALRecoUtils->SetEMCALChannelRecalibrationFactor(iSM,iCol,iRow,c)                                ; }
	
  void          SetPHOSChannelRecalibrationFactor (Int_t imod, Int_t iCol, Int_t iRow, Double_t c = 1) {
                  if(!fPHOSRecalibrationFactors)  InitPHOSRecalibrationFactors();
                  ((TH2F*)fPHOSRecalibrationFactors->At(imod))->SetBinContent(iCol,iRow,c)                            ; }
    
  void          SetEMCALChannelRecalibrationFactors(Int_t iSM , TH2F* h) { fEMCALRecoUtils->SetEMCALChannelRecalibrationFactors(iSM,h)      ; }
  void          SetPHOSChannelRecalibrationFactors(Int_t imod , TH2F* h) { fPHOSRecalibrationFactors ->AddAt(h,imod)                        ; }
	
  TH2F *        GetEMCALChannelRecalibrationFactors(Int_t iSM)     const { return fEMCALRecoUtils->GetEMCALChannelRecalibrationFactors(iSM) ; }
  TH2F *        GetPHOSChannelRecalibrationFactors(Int_t imod)     const { return (TH2F*)fPHOSRecalibrationFactors->At(imod)                ; }
	
  void          SetEMCALChannelRecalibrationFactors(TObjArray *map)      { fEMCALRecoUtils->SetEMCALChannelRecalibrationFactors(map)        ; }
  void          SetPHOSChannelRecalibrationFactors (TObjArray *map)      { fPHOSRecalibrationFactors  = map;}

  void          RecalibrateCellTime     (Double_t & time, Int_t calo, Int_t absId, Int_t bunchCrossNumber) const ;
  void          RecalibrateCellAmplitude(Float_t  & amp,  Int_t calo, Int_t absId) const ;
  Float_t       RecalibrateClusterEnergy(AliVCluster* cluster, AliVCaloCells * cells);
  Float_t       RecalibrateClusterEnergyWeightCell(AliVCluster* cluster, AliVCaloCells * cells, Float_t energyOrg);

  // Run dependent energy calibrations (EMCAL)
  
  void          SwitchOffRunDepCorrection()                              {  fRunDependentCorrection = kFALSE  ; }
  void          SwitchOnRunDepCorrection()                               {  fRunDependentCorrection = kTRUE   ; }
  
  // Time Recalibration (EMCAL)
  
  Bool_t       IsTimeRecalibrationOn()                             const { return fEMCALRecoUtils->IsTimeRecalibrationOn() ; }
  void         SwitchOffTimeRecalibration()                              { fEMCALRecoUtils->SwitchOffTimeRecalibration()   ; }
  void         SwitchOnTimeRecalibration()                               { fEMCALRecoUtils->SwitchOnTimeRecalibration()    ; }
  
  Float_t      GetEMCALChannelTimeRecalibrationFactor(Int_t bc, Int_t absID) const
  { return fEMCALRecoUtils->GetEMCALChannelTimeRecalibrationFactor(bc, absID) ; } 
	
  void         SetEMCALChannelTimeRecalibrationFactor(Int_t bc, Int_t absID, Double_t c = 0)
  { fEMCALRecoUtils->SetEMCALChannelTimeRecalibrationFactor(bc, absID, c) ; }  
  
  TH1F *       GetEMCALChannelTimeRecalibrationFactors(Int_t bc) const     { return fEMCALRecoUtils-> GetEMCALChannelTimeRecalibrationFactors(bc) ; }
  void         SetEMCALChannelTimeRecalibrationFactors(TObjArray *map)     { fEMCALRecoUtils->SetEMCALChannelTimeRecalibrationFactors(map)        ; }
  void         SetEMCALChannelTimeRecalibrationFactors(Int_t bc , TH1F* h) { fEMCALRecoUtils->SetEMCALChannelTimeRecalibrationFactors(bc , h)     ; }
  
  //EMCAL specific utils for the moment
  void          SetEMCALRecoUtils(AliEMCALRecoUtils * ru)  { fEMCALRecoUtils = ru          ; }
  AliEMCALRecoUtils* GetEMCALRecoUtils()             const { return fEMCALRecoUtils        ; }
  
  Bool_t        IsCorrectionOfClusterEnergyOn()      const { return fCorrectELinearity     ; }
  void          SwitchOnCorrectClusterLinearity()          { fCorrectELinearity = kTRUE    ; } 
  void          SwitchOffCorrectClusterLinearity()         { fCorrectELinearity = kFALSE   ; } 
  void          CorrectClusterEnergy(AliVCluster *cl);
  
  Bool_t        IsRecalculationOfClusterPositionOn() const { return fRecalculatePosition   ; }
  void          SwitchOnRecalculateClusterPosition()       { fRecalculatePosition = kTRUE  ; } 
  void          SwitchOffRecalculateClusterPosition()      { fRecalculatePosition = kFALSE ; } 
  void          RecalculateClusterPosition(AliVCaloCells* cells, AliVCluster* clu);
  void          RecalculateClusterShowerShapeParameters(AliVCaloCells* cells, AliVCluster* clu){
                  fEMCALRecoUtils->RecalculateClusterShowerShapeParameters((AliEMCALGeometry*)fEMCALGeo, cells, clu)  ; }
  
  void          RecalculateClusterDistanceToBadChannel(AliVCaloCells* cells, AliVCluster* clu){
                  fEMCALRecoUtils->RecalculateClusterDistanceToBadChannel((AliEMCALGeometry*)fEMCALGeo, cells, clu)   ; }
  
  void          RecalculateClusterPID(AliVCluster* clu)    { fEMCALRecoUtils->RecalculateClusterPID(clu)              ; }

  // *** Track Matching ***
  
  AliVTrack *   GetMatchedTrack(AliVCluster * cluster, AliVEvent * event, Int_t index = 0) const ;
  
  // Recalculation
  void          RecalculateClusterTrackMatching(AliVEvent * event, TObjArray* clusterArray = 0x0) ;
    
  void          GetMatchedResiduals(Int_t index, Float_t &dR, Float_t &dZ) {
                  if (fRecalculateMatching) fEMCALRecoUtils->GetMatchedResiduals(index,dR,dZ)   ; }
  
  //This could be used for PHOS ...
  void          SwitchOnRecalculateClusterTrackMatching()       { fRecalculateMatching = kTRUE  ; } 
  void          SwitchOffRecalculateClusterTrackMatching()      { fRecalculateMatching = kFALSE ; } 
  Bool_t        IsRecalculationOfClusterTrackMatchingOn() const { return fRecalculateMatching   ; }
  
  Float_t       GetCutZ()                                 const { return fCutZ                  ; } // PHOS only
  void          SetCutZ(Float_t z)                              { fCutZ = z                     ; } // PHOS only

  
  Float_t       GetCutR()                                 const { return fCutR                  ; } // PHOS and EMCAL
  void          SetCutR(Float_t r)                              { fCutR = r                     ;   // PHOS and EMCA
                                                                  fEMCALRecoUtils->SetCutR(r)   ; }
  
  Float_t       GetCutEta()                               const { return fCutEta                ; } // EMCAL only
  void          SetCutEta(Float_t e)                            { fCutEta = e                   ;   // EMCAL only
                                                                  fEMCALRecoUtils->SetCutEta(e) ; }

  Float_t       GetCutPhi()                               const { return fCutPhi                ; } // EMCAL only
  void          SetCutPhi(Float_t p)                            { fCutPhi = p                   ;   // EMCAL only
                                                                  fEMCALRecoUtils->SetCutPhi(p) ; }
  // OADB options settings
  
  void          AccessOADB(AliVEvent * event) ;
  
  TString       GetPass() ;
  
  void          SwitchOnEMCALOADB()                             { fOADBForEMCAL = kTRUE         ; }
  void          SwitchOffEMCALOADB()                            { fOADBForEMCAL = kFALSE        ; }

  void          SwitchOnPHOSOADB()                              { fOADBForPHOS  = kTRUE         ; }
  void          SwitchOffPHOSOADB()                             { fOADBForPHOS  = kFALSE        ; }

  void          SetEMCALOADBFilePath(TString path)              { fOADBFilePathEMCAL  = path    ; }
  void          SetPHOSOADBFilePath (TString path)              { fOADBFilePathPHOS   = path    ; }

  void          SetNumberOfSuperModulesUsed(Int_t nSM)          { fNSuperModulesUsed  = nSM     ; }
  Int_t         GetNumberOfSuperModulesUsed()             const { return fNSuperModulesUsed     ; }
  
  enum detector { kEMCAL = 0, kPHOS = 1, kCTS = 2, kDCAL = 3, kDCALPHOS = 4 };
  
 private:

  Int_t              fDebug;                 //  Debugging level
  TString            fEMCALGeoName;          //  Name of geometry to use for EMCAL.
  TString            fPHOSGeoName;           //  Name of geometry to use for PHOS.	
  AliEMCALGeometry * fEMCALGeo ;             //! EMCAL geometry pointer
  AliPHOSGeoUtils  * fPHOSGeo  ;             //! PHOS  geometry pointer  
  Bool_t             fEMCALGeoMatrixSet;     //  Check if the transformation matrix is set for EMCAL
  Bool_t             fPHOSGeoMatrixSet ;     //  Check if the transformation matrix is set for PHOS
  Bool_t             fLoadEMCALMatrices;     //  Matrices set from configuration, not get from geometry.root or from ESDs/AODs
  TGeoHMatrix *      fEMCALMatrix[12];       //  Geometry matrices with alignments
  Bool_t             fLoadPHOSMatrices;      //  Matrices set from configuration, not get from geometry.root or from ESDs/AODs
  TGeoHMatrix *      fPHOSMatrix[5];         //  Geometry matrices with alignments
  Bool_t             fRemoveBadChannels;     //  Check the channel status provided and remove clusters with bad channels
  TObjArray        * fPHOSBadChannelMap;     //  Array of histograms with map of bad channels, PHOS
  Int_t              fNCellsFromPHOSBorder;  //  Number of cells from PHOS  border the cell with maximum amplitude has to be.
  Int_t              fNMaskCellColumns;      //  Number of masked columns
  Int_t*             fMaskCellColumns;       //[fNMaskCellColumns] list of masked cell collumn
  Bool_t             fRecalibration;         //  Switch on or off the recalibration
  Bool_t             fRunDependentCorrection;//  Switch on or off the recalibration dependent on T
  TObjArray        * fPHOSRecalibrationFactors;  // Array of histograms with map of recalibration factors, PHOS
  AliEMCALRecoUtils* fEMCALRecoUtils;        //  EMCAL utils for cluster rereconstruction
  Bool_t             fRecalculatePosition;   //  Recalculate cluster position
  Bool_t             fCorrectELinearity  ;   //  Correct cluster energy linearity
  Bool_t             fRecalculateMatching;   //  Recalculate cluster position
  Float_t            fCutR;                  //  dR cut on matching (PHOS)
  Float_t            fCutZ;                  //  dZ cut on matching (EMCAL/PHOS)
  Float_t            fCutEta;                //  dEta cut on matching (EMCAL)
  Float_t            fCutPhi;                //  dPhi cut on matching (EMCAL)
  Float_t            fLocMaxCutE;            //  Local maxima cut must have more than this energy
  Float_t            fLocMaxCutEDiff;        //  Local maxima cut, when aggregating cells, next can be a bit higher
  Bool_t             fPlotCluster;           //  Plot cluster in splitting method
  Bool_t             fOADBSet ;              //  AODB parameters already set
  Bool_t             fOADBForEMCAL ;         //  Get calibration from OADB for EMCAL
  Bool_t             fOADBForPHOS ;          //  Get calibration from OADB for PHOS
  TString            fOADBFilePathEMCAL ;    //  Default path $ALICE_ROOT/OADB/EMCAL, if needed change
  TString            fOADBFilePathPHOS ;     //  Default path $ALICE_ROOT/OADB/PHOS, if needed change
  Bool_t             fImportGeometryFromFile;// Import geometry settings in geometry.root file
  TString            fImportGeometryFilePath;// path fo geometry.root file

  Int_t              fNSuperModulesUsed;     // Number of supermodules to be used in analysis, can be different than the real geo,
                                             // to be used at initialization of histograms

  Bool_t             fMCECellClusFracCorrOn; // Correct or not the weight of cells in cluster
  Float_t            fMCECellClusFracCorrParam[4]; // Parameters for the function correcting the weight of the cells in the cluster
  
  AliCalorimeterUtils(              const AliCalorimeterUtils & cu) ; // cpy ctor
  AliCalorimeterUtils & operator = (const AliCalorimeterUtils & cu) ; // cpy assignment
  
  ClassDef(AliCalorimeterUtils,17)
} ;


#endif //ALICALORIMETERUTILS_H



 AliCalorimeterUtils.h:1
 AliCalorimeterUtils.h:2
 AliCalorimeterUtils.h:3
 AliCalorimeterUtils.h:4
 AliCalorimeterUtils.h:5
 AliCalorimeterUtils.h:6
 AliCalorimeterUtils.h:7
 AliCalorimeterUtils.h:8
 AliCalorimeterUtils.h:9
 AliCalorimeterUtils.h:10
 AliCalorimeterUtils.h:11
 AliCalorimeterUtils.h:12
 AliCalorimeterUtils.h:13
 AliCalorimeterUtils.h:14
 AliCalorimeterUtils.h:15
 AliCalorimeterUtils.h:16
 AliCalorimeterUtils.h:17
 AliCalorimeterUtils.h:18
 AliCalorimeterUtils.h:19
 AliCalorimeterUtils.h:20
 AliCalorimeterUtils.h:21
 AliCalorimeterUtils.h:22
 AliCalorimeterUtils.h:23
 AliCalorimeterUtils.h:24
 AliCalorimeterUtils.h:25
 AliCalorimeterUtils.h:26
 AliCalorimeterUtils.h:27
 AliCalorimeterUtils.h:28
 AliCalorimeterUtils.h:29
 AliCalorimeterUtils.h:30
 AliCalorimeterUtils.h:31
 AliCalorimeterUtils.h:32
 AliCalorimeterUtils.h:33
 AliCalorimeterUtils.h:34
 AliCalorimeterUtils.h:35
 AliCalorimeterUtils.h:36
 AliCalorimeterUtils.h:37
 AliCalorimeterUtils.h:38
 AliCalorimeterUtils.h:39
 AliCalorimeterUtils.h:40
 AliCalorimeterUtils.h:41
 AliCalorimeterUtils.h:42
 AliCalorimeterUtils.h:43
 AliCalorimeterUtils.h:44
 AliCalorimeterUtils.h:45
 AliCalorimeterUtils.h:46
 AliCalorimeterUtils.h:47
 AliCalorimeterUtils.h:48
 AliCalorimeterUtils.h:49
 AliCalorimeterUtils.h:50
 AliCalorimeterUtils.h:51
 AliCalorimeterUtils.h:52
 AliCalorimeterUtils.h:53
 AliCalorimeterUtils.h:54
 AliCalorimeterUtils.h:55
 AliCalorimeterUtils.h:56
 AliCalorimeterUtils.h:57
 AliCalorimeterUtils.h:58
 AliCalorimeterUtils.h:59
 AliCalorimeterUtils.h:60
 AliCalorimeterUtils.h:61
 AliCalorimeterUtils.h:62
 AliCalorimeterUtils.h:63
 AliCalorimeterUtils.h:64
 AliCalorimeterUtils.h:65
 AliCalorimeterUtils.h:66
 AliCalorimeterUtils.h:67
 AliCalorimeterUtils.h:68
 AliCalorimeterUtils.h:69
 AliCalorimeterUtils.h:70
 AliCalorimeterUtils.h:71
 AliCalorimeterUtils.h:72
 AliCalorimeterUtils.h:73
 AliCalorimeterUtils.h:74
 AliCalorimeterUtils.h:75
 AliCalorimeterUtils.h:76
 AliCalorimeterUtils.h:77
 AliCalorimeterUtils.h:78
 AliCalorimeterUtils.h:79
 AliCalorimeterUtils.h:80
 AliCalorimeterUtils.h:81
 AliCalorimeterUtils.h:82
 AliCalorimeterUtils.h:83
 AliCalorimeterUtils.h:84
 AliCalorimeterUtils.h:85
 AliCalorimeterUtils.h:86
 AliCalorimeterUtils.h:87
 AliCalorimeterUtils.h:88
 AliCalorimeterUtils.h:89
 AliCalorimeterUtils.h:90
 AliCalorimeterUtils.h:91
 AliCalorimeterUtils.h:92
 AliCalorimeterUtils.h:93
 AliCalorimeterUtils.h:94
 AliCalorimeterUtils.h:95
 AliCalorimeterUtils.h:96
 AliCalorimeterUtils.h:97
 AliCalorimeterUtils.h:98
 AliCalorimeterUtils.h:99
 AliCalorimeterUtils.h:100
 AliCalorimeterUtils.h:101
 AliCalorimeterUtils.h:102
 AliCalorimeterUtils.h:103
 AliCalorimeterUtils.h:104
 AliCalorimeterUtils.h:105
 AliCalorimeterUtils.h:106
 AliCalorimeterUtils.h:107
 AliCalorimeterUtils.h:108
 AliCalorimeterUtils.h:109
 AliCalorimeterUtils.h:110
 AliCalorimeterUtils.h:111
 AliCalorimeterUtils.h:112
 AliCalorimeterUtils.h:113
 AliCalorimeterUtils.h:114
 AliCalorimeterUtils.h:115
 AliCalorimeterUtils.h:116
 AliCalorimeterUtils.h:117
 AliCalorimeterUtils.h:118
 AliCalorimeterUtils.h:119
 AliCalorimeterUtils.h:120
 AliCalorimeterUtils.h:121
 AliCalorimeterUtils.h:122
 AliCalorimeterUtils.h:123
 AliCalorimeterUtils.h:124
 AliCalorimeterUtils.h:125
 AliCalorimeterUtils.h:126
 AliCalorimeterUtils.h:127
 AliCalorimeterUtils.h:128
 AliCalorimeterUtils.h:129
 AliCalorimeterUtils.h:130
 AliCalorimeterUtils.h:131
 AliCalorimeterUtils.h:132
 AliCalorimeterUtils.h:133
 AliCalorimeterUtils.h:134
 AliCalorimeterUtils.h:135
 AliCalorimeterUtils.h:136
 AliCalorimeterUtils.h:137
 AliCalorimeterUtils.h:138
 AliCalorimeterUtils.h:139
 AliCalorimeterUtils.h:140
 AliCalorimeterUtils.h:141
 AliCalorimeterUtils.h:142
 AliCalorimeterUtils.h:143
 AliCalorimeterUtils.h:144
 AliCalorimeterUtils.h:145
 AliCalorimeterUtils.h:146
 AliCalorimeterUtils.h:147
 AliCalorimeterUtils.h:148
 AliCalorimeterUtils.h:149
 AliCalorimeterUtils.h:150
 AliCalorimeterUtils.h:151
 AliCalorimeterUtils.h:152
 AliCalorimeterUtils.h:153
 AliCalorimeterUtils.h:154
 AliCalorimeterUtils.h:155
 AliCalorimeterUtils.h:156
 AliCalorimeterUtils.h:157
 AliCalorimeterUtils.h:158
 AliCalorimeterUtils.h:159
 AliCalorimeterUtils.h:160
 AliCalorimeterUtils.h:161
 AliCalorimeterUtils.h:162
 AliCalorimeterUtils.h:163
 AliCalorimeterUtils.h:164
 AliCalorimeterUtils.h:165
 AliCalorimeterUtils.h:166
 AliCalorimeterUtils.h:167
 AliCalorimeterUtils.h:168
 AliCalorimeterUtils.h:169
 AliCalorimeterUtils.h:170
 AliCalorimeterUtils.h:171
 AliCalorimeterUtils.h:172
 AliCalorimeterUtils.h:173
 AliCalorimeterUtils.h:174
 AliCalorimeterUtils.h:175
 AliCalorimeterUtils.h:176
 AliCalorimeterUtils.h:177
 AliCalorimeterUtils.h:178
 AliCalorimeterUtils.h:179
 AliCalorimeterUtils.h:180
 AliCalorimeterUtils.h:181
 AliCalorimeterUtils.h:182
 AliCalorimeterUtils.h:183
 AliCalorimeterUtils.h:184
 AliCalorimeterUtils.h:185
 AliCalorimeterUtils.h:186
 AliCalorimeterUtils.h:187
 AliCalorimeterUtils.h:188
 AliCalorimeterUtils.h:189
 AliCalorimeterUtils.h:190
 AliCalorimeterUtils.h:191
 AliCalorimeterUtils.h:192
 AliCalorimeterUtils.h:193
 AliCalorimeterUtils.h:194
 AliCalorimeterUtils.h:195
 AliCalorimeterUtils.h:196
 AliCalorimeterUtils.h:197
 AliCalorimeterUtils.h:198
 AliCalorimeterUtils.h:199
 AliCalorimeterUtils.h:200
 AliCalorimeterUtils.h:201
 AliCalorimeterUtils.h:202
 AliCalorimeterUtils.h:203
 AliCalorimeterUtils.h:204
 AliCalorimeterUtils.h:205
 AliCalorimeterUtils.h:206
 AliCalorimeterUtils.h:207
 AliCalorimeterUtils.h:208
 AliCalorimeterUtils.h:209
 AliCalorimeterUtils.h:210
 AliCalorimeterUtils.h:211
 AliCalorimeterUtils.h:212
 AliCalorimeterUtils.h:213
 AliCalorimeterUtils.h:214
 AliCalorimeterUtils.h:215
 AliCalorimeterUtils.h:216
 AliCalorimeterUtils.h:217
 AliCalorimeterUtils.h:218
 AliCalorimeterUtils.h:219
 AliCalorimeterUtils.h:220
 AliCalorimeterUtils.h:221
 AliCalorimeterUtils.h:222
 AliCalorimeterUtils.h:223
 AliCalorimeterUtils.h:224
 AliCalorimeterUtils.h:225
 AliCalorimeterUtils.h:226
 AliCalorimeterUtils.h:227
 AliCalorimeterUtils.h:228
 AliCalorimeterUtils.h:229
 AliCalorimeterUtils.h:230
 AliCalorimeterUtils.h:231
 AliCalorimeterUtils.h:232
 AliCalorimeterUtils.h:233
 AliCalorimeterUtils.h:234
 AliCalorimeterUtils.h:235
 AliCalorimeterUtils.h:236
 AliCalorimeterUtils.h:237
 AliCalorimeterUtils.h:238
 AliCalorimeterUtils.h:239
 AliCalorimeterUtils.h:240
 AliCalorimeterUtils.h:241
 AliCalorimeterUtils.h:242
 AliCalorimeterUtils.h:243
 AliCalorimeterUtils.h:244
 AliCalorimeterUtils.h:245
 AliCalorimeterUtils.h:246
 AliCalorimeterUtils.h:247
 AliCalorimeterUtils.h:248
 AliCalorimeterUtils.h:249
 AliCalorimeterUtils.h:250
 AliCalorimeterUtils.h:251
 AliCalorimeterUtils.h:252
 AliCalorimeterUtils.h:253
 AliCalorimeterUtils.h:254
 AliCalorimeterUtils.h:255
 AliCalorimeterUtils.h:256
 AliCalorimeterUtils.h:257
 AliCalorimeterUtils.h:258
 AliCalorimeterUtils.h:259
 AliCalorimeterUtils.h:260
 AliCalorimeterUtils.h:261
 AliCalorimeterUtils.h:262
 AliCalorimeterUtils.h:263
 AliCalorimeterUtils.h:264
 AliCalorimeterUtils.h:265
 AliCalorimeterUtils.h:266
 AliCalorimeterUtils.h:267
 AliCalorimeterUtils.h:268
 AliCalorimeterUtils.h:269
 AliCalorimeterUtils.h:270
 AliCalorimeterUtils.h:271
 AliCalorimeterUtils.h:272
 AliCalorimeterUtils.h:273
 AliCalorimeterUtils.h:274
 AliCalorimeterUtils.h:275
 AliCalorimeterUtils.h:276
 AliCalorimeterUtils.h:277
 AliCalorimeterUtils.h:278
 AliCalorimeterUtils.h:279
 AliCalorimeterUtils.h:280
 AliCalorimeterUtils.h:281
 AliCalorimeterUtils.h:282
 AliCalorimeterUtils.h:283
 AliCalorimeterUtils.h:284
 AliCalorimeterUtils.h:285
 AliCalorimeterUtils.h:286
 AliCalorimeterUtils.h:287
 AliCalorimeterUtils.h:288
 AliCalorimeterUtils.h:289
 AliCalorimeterUtils.h:290
 AliCalorimeterUtils.h:291
 AliCalorimeterUtils.h:292
 AliCalorimeterUtils.h:293
 AliCalorimeterUtils.h:294
 AliCalorimeterUtils.h:295
 AliCalorimeterUtils.h:296
 AliCalorimeterUtils.h:297
 AliCalorimeterUtils.h:298
 AliCalorimeterUtils.h:299
 AliCalorimeterUtils.h:300
 AliCalorimeterUtils.h:301
 AliCalorimeterUtils.h:302
 AliCalorimeterUtils.h:303
 AliCalorimeterUtils.h:304
 AliCalorimeterUtils.h:305
 AliCalorimeterUtils.h:306
 AliCalorimeterUtils.h:307
 AliCalorimeterUtils.h:308
 AliCalorimeterUtils.h:309
 AliCalorimeterUtils.h:310
 AliCalorimeterUtils.h:311
 AliCalorimeterUtils.h:312
 AliCalorimeterUtils.h:313
 AliCalorimeterUtils.h:314
 AliCalorimeterUtils.h:315
 AliCalorimeterUtils.h:316
 AliCalorimeterUtils.h:317
 AliCalorimeterUtils.h:318
 AliCalorimeterUtils.h:319
 AliCalorimeterUtils.h:320
 AliCalorimeterUtils.h:321
 AliCalorimeterUtils.h:322
 AliCalorimeterUtils.h:323
 AliCalorimeterUtils.h:324
 AliCalorimeterUtils.h:325
 AliCalorimeterUtils.h:326
 AliCalorimeterUtils.h:327
 AliCalorimeterUtils.h:328
 AliCalorimeterUtils.h:329
 AliCalorimeterUtils.h:330
 AliCalorimeterUtils.h:331
 AliCalorimeterUtils.h:332
 AliCalorimeterUtils.h:333
 AliCalorimeterUtils.h:334
 AliCalorimeterUtils.h:335
 AliCalorimeterUtils.h:336
 AliCalorimeterUtils.h:337
 AliCalorimeterUtils.h:338
 AliCalorimeterUtils.h:339
 AliCalorimeterUtils.h:340
 AliCalorimeterUtils.h:341
 AliCalorimeterUtils.h:342
 AliCalorimeterUtils.h:343
 AliCalorimeterUtils.h:344
 AliCalorimeterUtils.h:345
 AliCalorimeterUtils.h:346
 AliCalorimeterUtils.h:347
 AliCalorimeterUtils.h:348
 AliCalorimeterUtils.h:349
 AliCalorimeterUtils.h:350
 AliCalorimeterUtils.h:351
 AliCalorimeterUtils.h:352
 AliCalorimeterUtils.h:353
 AliCalorimeterUtils.h:354
 AliCalorimeterUtils.h:355
 AliCalorimeterUtils.h:356
 AliCalorimeterUtils.h:357
 AliCalorimeterUtils.h:358
 AliCalorimeterUtils.h:359
 AliCalorimeterUtils.h:360
 AliCalorimeterUtils.h:361
 AliCalorimeterUtils.h:362
 AliCalorimeterUtils.h:363
 AliCalorimeterUtils.h:364
 AliCalorimeterUtils.h:365
 AliCalorimeterUtils.h:366
 AliCalorimeterUtils.h:367
 AliCalorimeterUtils.h:368
 AliCalorimeterUtils.h:369
 AliCalorimeterUtils.h:370
 AliCalorimeterUtils.h:371