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

//_________________________________________________________________________
// Class for PID selection with calorimeters
// The Output of the main method GetIdentifiedParticleType is a PDG number identifying the cluster, 
// being kPhoton, kElectron, kPi0 ... as defined in the header file
//   - GetIdentifiedParticleType(const AliVCluster * cluster) 
//      Assignes a PID tag to the cluster, right now there is the possibility to : use bayesian weights from reco, 
//      recalculate them (EMCAL) or use other procedures not used in reco.
//      In order to recalculate Bayesian, it is necessary to load the EMCALUtils library
//      and do SwitchOnBayesianRecalculation().
//      To change the PID parameters from Low to High like the ones by default, use the constructor 
//      AliCaloPID(flux)
//      where flux is AliCaloPID::kLow or AliCaloPID::kHigh
//      If it is necessary to change the parameters use the constructor 
//      AliCaloPID(AliEMCALPIDUtils *utils) and set the parameters before.

//   - GetGetIdentifiedParticleTypeFromBayesian(const Double_t * pid, const Float_t energy)
//      Reads the PID weights array of the ESDs and depending on its magnitude identifies the particle, 
//      executed when bayesian is ON by GetIdentifiedParticleType(const AliVCluster * cluster) 
//   - SetPIDBits: Simple PID, depending on the thresholds fLOCut fTOFCut and even the
//     result of the PID bayesian a different PID bit is set. 
//
//
//*-- Author: Gustavo Conesa (INFN-LNF)

// --- ROOT system ---
#include <TObject.h> 
class TString ;
class TLorentzVector ;
#include <TFormula.h>
class TList;
class TH2F ;

//--- AliRoot system ---
class AliVCluster;
class AliVCaloCells;
class AliAODPWG4Particle;
class AliEMCALPIDUtils;
class AliCalorimeterUtils;
class AliVEvent;

class AliCaloPID : public TObject {
	
 public: 
  
  AliCaloPID() ; // ctor
  AliCaloPID(Int_t particleFlux) ; // ctor, to be used when recalculating bayesian PID
  AliCaloPID(const TNamed * emcalpid) ; // ctor, to be used when recalculating bayesian PID and need different parameters
  virtual ~AliCaloPID() ;//virtual dtor
  	
  enum PidType 
  {
    kPhoton         = 22,
    kPi0            = 111,
    kEta            = 221, 
    kElectron       = 11, 
    kEleCon         =-11, 
    kNeutralHadron  = 2112, 
    kChargedHadron  = 211, 
    kNeutralUnknown = 130, 
    kChargedUnknown = 321
  };
  
  enum TagType {kPi0Decay, kEtaDecay, kOtherDecay, kConversion, kNoTag = -1};
  
  // Main methods
  
  TList *   GetCreateOutputObjects();

  void      InitParameters();
  
  Bool_t    IsInPi0SplitAsymmetryRange(Float_t energy, Float_t asy,  Int_t nlm) const;

  Bool_t    IsInPi0SplitMassRange     (Float_t energy, Float_t mass, Int_t nlm) const;
  
  Bool_t    IsInM02Range              (Float_t m02) const;
  Bool_t    IsInPi0M02Range           (Float_t energy, Float_t m02,  Int_t nlm) const;
  Bool_t    IsInEtaM02Range           (Float_t energy, Float_t m02,  Int_t nlm) const;
  Bool_t    IsInConM02Range           (Float_t energy, Float_t m02,  Int_t nlm) const;
  
  
  Int_t     GetIdentifiedParticleTypeFromBayesWeights(Bool_t isEMCAL, Double_t * pid, Float_t energy) ;

  Int_t     GetIdentifiedParticleTypeFromClusterSplitting(AliVCluster * cluster, AliVCaloCells* cells, 
                                                          AliCalorimeterUtils * caloutils,
                                                          Double_t vertex[3], 
                                                          Int_t & nLocMax, Double_t & mass, Double_t & angle,
                                                          TLorentzVector & l1  , TLorentzVector & l2,
                                                          Int_t   & absId1,   Int_t   & absId2,
                                                          Float_t & distbad1, Float_t & distbad2,
                                                          Bool_t  & fidcut1,  Bool_t  & fidcut2  ) const;
  
  Int_t     GetIdentifiedParticleType(AliVCluster * cluster) ;
  
  TString   GetPIDParametersList();
  
  Bool_t    IsTrackMatched(AliVCluster * cluster, AliCalorimeterUtils* cu, AliVEvent* event) const ;    
  
  void      SetPIDBits(AliVCluster * cluster, AliAODPWG4Particle *aodph, 
                       AliCalorimeterUtils* cu, AliVEvent* event);
  
  void      Print(const Option_t * opt)const;
  
  void      PrintClusterPIDWeights(const Double_t * pid) const;
  
  //Check if cluster photon-like. Uses photon cluster parameterization in real pp data 
  //Returns distance in sigmas. Recommended cut 2.5
  Float_t TestPHOSDispersion(Double_t pt, Double_t m20, Double_t m02) const ;
  //Checks distance to the closest track. Takes into account 
  //non-perpendicular incidence of tracks.
  Float_t TestPHOSChargedVeto(Double_t dx, Double_t dz, Double_t ptTrack,
                              Int_t chargeTrack, Double_t mf) const ;
  
  // Setters, getters
  
  void    SetDebug(Int_t deb)                  { fDebug = deb                 ; }
  Int_t   GetDebug()                     const { return fDebug                ; }	

  enum    eventType{kLow,kHigh};
  void    SetLowParticleFlux()                 { fParticleFlux        = kLow  ; }
  void    SetHighParticleFlux()                { fParticleFlux        = kHigh ; }  
  // not really used, only for bayesian recalculation in EMCAL, but could be useful in future
  
  // Bayesian
  
  void    SwitchOnBayesian()                   { fUseBayesianWeights  = kTRUE ; }
  void    SwitchOffBayesian()                  { fUseBayesianWeights  = kFALSE; }
  void    SwitchOnBayesianRecalculation()      { fRecalculateBayesian = kTRUE ; fUseBayesianWeights  = kTRUE ;} // EMCAL
  void    SwitchOffBayesianRecalculation()     { fRecalculateBayesian = kFALSE; } // EMCAL
  
  AliEMCALPIDUtils * GetEMCALPIDUtils() ; 
  
  //Weight getters
  Float_t GetEMCALPhotonWeight()         const { return fEMCALPhotonWeight    ; }
  Float_t GetEMCALPi0Weight()            const { return fEMCALPi0Weight       ; }
  Float_t GetEMCALElectronWeight()       const { return fEMCALElectronWeight  ; }
  Float_t GetEMCALChargeWeight()         const { return fEMCALChargeWeight    ; }
  Float_t GetEMCALNeutralWeight()        const { return fEMCALNeutralWeight   ; }
  Float_t GetPHOSPhotonWeight()          const { return fPHOSPhotonWeight     ; }
  Float_t GetPHOSPi0Weight()             const { return fPHOSPi0Weight        ; }
  Float_t GetPHOSElectronWeight()        const { return fPHOSElectronWeight   ; }
  Float_t GetPHOSChargeWeight()          const { return fPHOSChargeWeight     ; }
  Float_t GetPHOSNeutralWeight()         const { return fPHOSNeutralWeight    ; }
  
  Bool_t  IsPHOSPIDWeightFormulaOn()     const { return fPHOSWeightFormula    ; } 

  TFormula * GetPHOSPhotonWeightFormula()      { 
    if(!fPHOSPhotonWeightFormula) 
      fPHOSPhotonWeightFormula = new TFormula("phos_photon_weight",
                                              fPHOSPhotonWeightFormulaExpression);
    return fPHOSPhotonWeightFormula                                           ; } 
  
  TFormula * GetPHOSPi0WeightFormula()         { 
    if(!fPHOSPi0WeightFormula) 
      fPHOSPi0WeightFormula = new TFormula("phos_pi0_weight",
                                           fPHOSPi0WeightFormulaExpression);
    return fPHOSPi0WeightFormula                                              ; } 
  
  TString  GetPHOSPhotonWeightFormulaExpression() const { return fPHOSPhotonWeightFormulaExpression ; } 
  TString  GetPHOSPi0WeightFormulaExpression()    const { return fPHOSPi0WeightFormulaExpression    ; } 
  
  //Weight setters
  void    SetEMCALPhotonWeight  (Float_t  w)   { fEMCALPhotonWeight   = w     ; }
  void    SetEMCALPi0Weight     (Float_t  w)   { fEMCALPi0Weight      = w     ; }
  void    SetEMCALElectronWeight(Float_t  w)   { fEMCALElectronWeight = w     ; }
  void    SetEMCALChargeWeight  (Float_t  w)   { fEMCALChargeWeight   = w     ; }
  void    SetEMCALNeutralWeight (Float_t  w)   { fEMCALNeutralWeight  = w     ; }
  void    SetPHOSPhotonWeight   (Float_t  w)   { fPHOSPhotonWeight    = w     ; }
  void    SetPHOSPi0Weight      (Float_t  w)   { fPHOSPi0Weight       = w     ; }
  void    SetPHOSElectronWeight (Float_t  w)   { fPHOSElectronWeight  = w     ; }
  void    SetPHOSChargeWeight   (Float_t  w)   { fPHOSChargeWeight    = w     ; }
  void    SetPHOSNeutralWeight  (Float_t  w)   { fPHOSNeutralWeight   = w     ; }
  
  void    UsePHOSPIDWeightFormula   (Bool_t ok )           { fPHOSWeightFormula                 = ok ; } 
  void    SetPHOSPhotonWeightFormulaExpression(TString ph) { fPHOSPhotonWeightFormulaExpression = ph ; } 
  void    SetPHOSPi0WeightFormulaExpression   (TString pi) { fPHOSPi0WeightFormulaExpression    = pi ; }
  
  //PID cuts 
  
  void    SetEMCALLambda0CutMax(Float_t lcut ) { fEMCALL0CutMax     = lcut    ; }
  Float_t GetEMCALLambda0CutMax()        const { return fEMCALL0CutMax        ; }   
  
  void    SetEMCALLambda0CutMin(Float_t lcut ) { fEMCALL0CutMin     = lcut    ; }
  Float_t GetEMCALLambda0CutMin()        const { return fEMCALL0CutMin        ; }   
  
  void    SetEMCALDEtaCut(Float_t dcut )       { fEMCALDEtaCut      = dcut    ; }
  Float_t GetEMCALDEtaCut()              const { return fEMCALDEtaCut         ; }   
  
  void    SetEMCALDPhiCut(Float_t dcut )       { fEMCALDPhiCut      = dcut    ; }
  Float_t GetEMCALDPhiCut()              const { return fEMCALDPhiCut         ; }   
  
  void    SetTOFCut(Float_t tcut )             { fTOFCut            = tcut    ; }
  Float_t GetTOFCut()                    const { return fTOFCut               ; }   
  
  void    SetPHOSRCut(Float_t rcut )           { fPHOSRCut          = rcut    ; }
  Float_t GetPHOSRCut()                  const { return fPHOSRCut             ; }   

  void    SetPHOSDispersionCut(Float_t dcut )  { fPHOSDispersionCut = dcut    ; }
  Float_t GetPHOSDispersionCut()         const { return fPHOSDispersionCut    ; }   
  
  // Cluster splitting analysis
  
  void    SwitchOnSimpleSplitMassCut()         { fUseSimpleMassCut   = kTRUE  ; }
  void    SwitchOffSimpleSplitMassCut()        { fUseSimpleMassCut   = kFALSE ; }

  void    SwitchOnSimpleSplitM02Cut()          { fUseSimpleM02Cut    = kTRUE  ; }
  void    SwitchOffSimpleSplitM02Cut()         { fUseSimpleM02Cut    = kFALSE ; }
  
  void    SwitchOnSplitAsymmetryCut()          { fUseSplitAsyCut     = kTRUE  ; }
  void    SwitchOffSplitAsymmetryCut()         { fUseSplitAsyCut     = kFALSE ; }
  Bool_t  IsSplitAsymmetryCutOn()              { return fUseSplitAsyCut       ; }

  void    SwitchOnSplitShowerShapeCut()        { fUseSplitSSCut      = kTRUE  ; }
  void    SwitchOffSplitShowerShapeCut()       { fUseSplitSSCut      = kFALSE ; }
  Bool_t  IsSplitShowerShapeCutOn()            { return fUseSplitSSCut        ; }
  
  void    SetClusterSplittingM02Cut(Float_t min=0, Float_t max=100) 
  { fSplitM02MinCut   = min ; fSplitM02MaxCut  = max ; }
  
  void    SetClusterSplittingMinNCells(Int_t c) { fSplitMinNCells = c         ; }
  Int_t   GetClusterSplittingMinNCells() const  { return fSplitMinNCells      ; }
  
  void    SetSplitEnergyFractionMinimum(Int_t i, Float_t min){ if (i < 3 && i >=0 ) fSplitEFracMin[i]  = min ; }
  Float_t GetSplitEnergyFractionMinimum(Int_t i) const       { if( i < 3 && i >=0 ) return fSplitEFracMin[i] ;  else return 0 ; }

  void    SetSubClusterEnergyMinimum   (Int_t i, Float_t min){ if (i < 3 && i >=0 ) fSubClusterEMin[i] = min ; }
  Float_t GetSubClusterEnergyMinimum   (Int_t i) const       { if( i < 3 && i >=0 ) return fSubClusterEMin[i];  else return 0 ; }
  
  Float_t GetPi0MinMass()                const { return fMassPi0Min           ; } // Simple cut case
  Float_t GetEtaMinMass()                const { return fMassEtaMin           ; } // Simple cut case
  Float_t GetPhotonMinMass()             const { return fMassPhoMin           ; }  
  Float_t GetPi0MaxMass()                const { return fMassPi0Max           ; }
  Float_t GetEtaMaxMass()                const { return fMassEtaMax           ; }
  Float_t GetPhotonMaxMass()             const { return fMassPhoMax           ; }
  
  void    SetSplitWidthSigma(Float_t s)        { fSplitWidthSigma        = s  ; }

  void    SetPi0MassShiftHighECell(Float_t s)  { fMassShiftHighECell     = s  ; }
  
  void    SetPi0MassSelectionParameters    (Int_t inlm, Int_t iparam, Float_t param)
  { if(iparam < 6 ) fMassPi0Param[inlm][iparam] = param ; }

  void    SetPi0WidthSelectionParameters    (Int_t inlm, Int_t iparam, Float_t param)
  { if(iparam < 6 ) fWidthPi0Param[inlm][iparam] = param ; }
  
  void    SetM02MaximumSelectionParameters      (Int_t inlm, Int_t iparam, Float_t param)
  { if(iparam < 5 && inlm < 2) fM02MaxParam[inlm][iparam] = param ; }
  
  void    SetM02MaximumShiftForNLMN(Int_t shift) { fM02MaxParamShiftNLMN = shift ; }
  
  void    SetM02MinimumSelectionParameters      (Int_t inlm, Int_t iparam, Float_t param)
  { if(iparam < 5 && inlm < 2) fM02MinParam[inlm][iparam] = param ; }
  
  void    SetAsymmetryMinimumSelectionParameters(Int_t inlm, Int_t iparam, Float_t param)
  { if(iparam < 4 && inlm < 2) fAsyMinParam[inlm][iparam] = param ; }

  void    SetPi0MassRange(Float_t min, Float_t max)    { fMassPi0Min    = min ; fMassPi0Max = max ; } // Simple case
  void    SetEtaMassRange(Float_t min, Float_t max)    { fMassEtaMin    = min ; fMassEtaMax = max ; }
  void    SetPhotonMassRange(Float_t min, Float_t max) { fMassPhoMin    = min ; fMassPhoMax = max ; }
    
private:
  
  Int_t	    fDebug;                             // Debug level
  Int_t     fParticleFlux;                      // Particle flux for setting PID parameters

  // Bayesian
  AliEMCALPIDUtils * fEMCALPIDUtils;            // Pointer to EMCALPID to redo the PID Bayesian calculation
  Bool_t    fUseBayesianWeights;                // Select clusters based on weights calculated in reconstruction
  Bool_t    fRecalculateBayesian;               // Recalculate PID bayesian or use simple PID?

  Float_t   fEMCALPhotonWeight;                 // Bayesian PID weight for photons in EMCAL 
  Float_t   fEMCALPi0Weight;                    // Bayesian PID weight for pi0 in EMCAL 
  Float_t   fEMCALElectronWeight;               // Bayesian PID weight for electrons in EMCAL 
  Float_t   fEMCALChargeWeight;                 // Bayesian PID weight for charged hadrons in EMCAL 
  Float_t   fEMCALNeutralWeight;                // Bayesian PID weight for neutral hadrons in EMCAL 
  Float_t   fPHOSPhotonWeight;                  // Bayesian PID weight for photons in PHOS 
  Float_t   fPHOSPi0Weight;                     // Bayesian PID weight for pi0 in PHOS 
  Float_t   fPHOSElectronWeight;                // Bayesian PID weight for electrons in PHOS 
  Float_t   fPHOSChargeWeight;                  // Bayesian PID weight for charged hadrons in PHOS 
  Float_t   fPHOSNeutralWeight;                 // Bayesian PID weight for neutral hadrons in PHOS 
  
  Bool_t    fPHOSWeightFormula ;                // Use parametrized weight threshold, function of energy
  TFormula *fPHOSPhotonWeightFormula ;          // Formula for photon weight
  TFormula *fPHOSPi0WeightFormula ;             // Formula for pi0 weight
  TString   fPHOSPhotonWeightFormulaExpression; // Photon weight formula in string
  TString   fPHOSPi0WeightFormulaExpression;    // Pi0 weight formula in string

  // PID calculation
  Float_t   fEMCALL0CutMax;                     // Max Cut on shower shape lambda0, used in PID evaluation, only EMCAL
  Float_t   fEMCALL0CutMin;                     // Min Cut on shower shape lambda0, used in PID evaluation, only EMCAL
  Float_t   fEMCALDEtaCut;                      // Track matching cut on Dz
  Float_t   fEMCALDPhiCut;                      // Track matching cut on Dx

  Float_t   fTOFCut;                            // Cut on TOF, used in PID evaluation
  
  Float_t   fPHOSDispersionCut;                 // Shower shape elipse radious cut
  Float_t   fPHOSRCut;                          // Track-Cluster distance cut for track matching in PHOS  
  
  // Cluster splitting mass ranges
  Bool_t    fUseSimpleMassCut;                  // Use simple min-max pi0 mass cut
  Bool_t    fUseSimpleM02Cut;                   // Use simple min-max M02 cut
  Bool_t    fUseSplitAsyCut ;                   // Remove splitted clusters with too large asymmetry
  Bool_t    fUseSplitSSCut  ;                   // Remove splitted clusters out of shower shape band
  Float_t   fSplitM02MaxCut ;                   // Study clusters with l0 smaller than cut
  Float_t   fSplitM02MinCut ;                   // Study clusters with l0 larger than cut  // simple case
  Int_t     fSplitMinNCells ;                   // Study clusters with ncells larger than cut  
  Float_t   fMassEtaMin  ;                      // Min Eta mass
  Float_t   fMassEtaMax  ;                      // Max Eta mass  
  Float_t   fMassPi0Min  ;                      // Min Pi0 mass // simple cut case
  Float_t   fMassPi0Max  ;                      // Min Pi0 mass // simple cut case
  Float_t   fMassPhoMin  ;                      // Min Photon mass
  Float_t   fMassPhoMax  ;                      // Min Photon mass
  Float_t   fMassPi0Param [2][6] ;              // mean mass param, 2 regions in energy
  Float_t   fWidthPi0Param[2][6] ;              // width param, 2 regions in energy
  Float_t   fM02MinParam[2][5] ;                // 5 param for expo + pol fit on M02 minimum for pi0 selection (maximum for conversions)
  Float_t   fM02MaxParam[2][5] ;                // 5 param for expo + pol fit on M02 maximum for pi0 selection
  Float_t   fM02MaxParamShiftNLMN;              // shift of max M02 for NLM>2
  Float_t   fAsyMinParam[2][4] ;                // 3 param for fit on asymmetry minimum, for 2 cases, NLM=1 and NLM>=2
  Float_t   fSplitEFracMin[3]  ;                // Do not use clusters with too large energy in cluster compared
                                                // to energy in splitted clusters, depeding on NLM
  Float_t   fSubClusterEMin[3]  ;               // Do not use sub-clusters with too low energy depeding on NLM
  Float_t   fSplitWidthSigma;                   // Cut on mass+-width*fSplitWidthSigma
  Float_t   fMassShiftHighECell;                // Shift cuts 5 MeV for Ecell > 150 MeV, default Ecell > 50 MeV

  AliCaloPID & operator = (const AliCaloPID & cpid) ; // cpy assignment
  AliCaloPID(              const AliCaloPID & cpid) ; // cpy ctor
  
  ClassDef(AliCaloPID,22)
  
} ;


#endif //ALICALOPID_H



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