]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/CaloTrackCorrBase/AliCaloPID.h
In case no maxima found because 2 high energy cells too close with energy difference...
[u/mrichter/AliRoot.git] / PWG / CaloTrackCorrBase / AliCaloPID.h
1 #ifndef ALICALOPID_H
2 #define ALICALOPID_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice     */
5
6 //_________________________________________________________________________
7 // Class for PID selection with calorimeters
8 // The Output of the main method GetIdentifiedParticleType is a PDG number identifying the cluster, 
9 // being kPhoton, kElectron, kPi0 ... as defined in the header file
10 //   - GetIdentifiedParticleType(const TString calo, const TLorentzVector mom, const AliVCluster * cluster) 
11 //      Assignes a PID tag to the cluster, right now there is the possibility to : use bayesian weights from reco, 
12 //      recalculate them (EMCAL) or use other procedures not used in reco.
13 //      In order to recalculate Bayesian, it is necessary to load the EMCALUtils library
14 //      and do SwitchOnBayesianRecalculation().
15 //      To change the PID parameters from Low to High like the ones by default, use the constructor 
16 //      AliCaloPID(flux)
17 //      where flux is AliCaloPID::kLow or AliCaloPID::kHigh
18 //      If it is necessary to change the parameters use the constructor 
19 //      AliCaloPID(AliEMCALPIDUtils *utils) and set the parameters before.
20
21 //   - GetGetIdentifiedParticleTypeFromBayesian(const TString calo, const Double_t * pid, const Float_t energy)
22 //      Reads the PID weights array of the ESDs and depending on its magnitude identifies the particle, 
23 //      executed when bayesian is ON by GetIdentifiedParticleType(const TString calo, const TLorentzVector mom, const AliVCluster * cluster) 
24 //   - SetPIDBits: Simple PID, depending on the thresholds fLOCut fTOFCut and even the
25 //     result of the PID bayesian a different PID bit is set. 
26 //
27 //
28 //*-- Author: Gustavo Conesa (INFN-LNF)
29
30 // --- ROOT system ---
31 #include <TObject.h> 
32 class TString ;
33 class TLorentzVector ;
34 #include <TFormula.h>
35 class TList;
36 class TH2F ;
37
38 //--- AliRoot system ---
39 class AliVCluster;
40 class AliAODPWG4Particle;
41 class AliEMCALPIDUtils;
42 class AliCalorimeterUtils;
43 class AliVEvent;
44
45 class AliCaloPID : public TObject {
46         
47  public: 
48   
49   AliCaloPID() ; // ctor
50   AliCaloPID(const Int_t particleFlux) ; // ctor, to be used when recalculating bayesian PID
51   AliCaloPID(const TNamed * emcalpid) ; // ctor, to be used when recalculating bayesian PID and need different parameters
52   virtual ~AliCaloPID() ;//virtual dtor
53         
54   enum PidType {
55     kPhoton         = 22,
56     kPi0            = 111,
57     kEta            = 221, 
58     kElectron       = 11, 
59     kEleCon         =-11, 
60     kNeutralHadron  = 2112, 
61     kChargedHadron  = 211, 
62     kNeutralUnknown = 130, 
63     kChargedUnknown = 321
64   };
65   
66   enum TagType {kPi0Decay, kEtaDecay, kOtherDecay, kConversion, kNoTag = -1};
67   
68   // Main methods
69   
70   TList *   GetCreateOutputObjects();
71
72   void      InitParameters();
73   
74   Int_t     GetIdentifiedParticleTypeFromBayesWeights(const TString calo, const Double_t * pid, const Float_t energy) ;
75   
76   Int_t     GetIdentifiedParticleType(const TString calo, const TLorentzVector mom, const AliVCluster * cluster) ;
77   
78   TString   GetPIDParametersList();
79   
80   Bool_t    IsTrackMatched(AliVCluster * cluster, AliCalorimeterUtils* cu, AliVEvent* event) const ;    
81   
82   void      SetPIDBits(const TString calo,  AliVCluster * cluster, AliAODPWG4Particle *aodph, 
83                        AliCalorimeterUtils* cu, AliVEvent* event);
84   
85   void      Print(const Option_t * opt)const;
86   
87   //Check if cluster photon-like. Uses photon cluster parameterization in real pp data 
88   //Returns distance in sigmas. Recommended cut 2.5
89   Float_t TestPHOSDispersion(const Double_t pt, const Double_t m20, const Double_t m02) const ; 
90   //Checks distance to the closest track. Takes into account 
91   //non-perpendicular incidence of tracks.
92   Float_t TestPHOSChargedVeto(const Double_t dx,  const Double_t dz, const Double_t ptTrack, 
93                               const Int_t chargeTrack, const Double_t mf) const ;
94   
95   // Setters, getters
96   
97   void    SetDebug(Int_t deb)                  { fDebug = deb                 ; }
98   Int_t   GetDebug()                     const { return fDebug                ; }       
99
100   enum    eventType{kLow,kHigh};
101   void    SetLowParticleFlux()                 { fParticleFlux        = kLow  ; }
102   void    SetHighParticleFlux()                { fParticleFlux        = kHigh ; }  
103   // not really used, only for bayesian recalculation in EMCAL, but could be useful in future
104   
105   // Bayesian
106   
107   void    SwitchOnBayesian()                   { fUseBayesianWeights  = kTRUE ; }
108   void    SwitchOffBayesian()                  { fUseBayesianWeights  = kFALSE; }
109   void    SwitchOnBayesianRecalculation()      { fRecalculateBayesian = kTRUE ; fUseBayesianWeights  = kTRUE ;} // EMCAL
110   void    SwitchOffBayesianRecalculation()     { fRecalculateBayesian = kFALSE; } // EMCAL
111   
112   AliEMCALPIDUtils * GetEMCALPIDUtils() ; 
113   
114   //Weight getters
115   Float_t GetEMCALPhotonWeight()         const { return fEMCALPhotonWeight    ; }
116   Float_t GetEMCALPi0Weight()            const { return fEMCALPi0Weight       ; }
117   Float_t GetEMCALElectronWeight()       const { return fEMCALElectronWeight  ; }
118   Float_t GetEMCALChargeWeight()         const { return fEMCALChargeWeight    ; }
119   Float_t GetEMCALNeutralWeight()        const { return fEMCALNeutralWeight   ; }
120   Float_t GetPHOSPhotonWeight()          const { return fPHOSPhotonWeight     ; }
121   Float_t GetPHOSPi0Weight()             const { return fPHOSPi0Weight        ; }
122   Float_t GetPHOSElectronWeight()        const { return fPHOSElectronWeight   ; }
123   Float_t GetPHOSChargeWeight()          const { return fPHOSChargeWeight     ; }
124   Float_t GetPHOSNeutralWeight()         const { return fPHOSNeutralWeight    ; }
125   
126   Bool_t  IsPHOSPIDWeightFormulaOn()     const { return fPHOSWeightFormula    ; } 
127
128   TFormula * GetPHOSPhotonWeightFormula()      { 
129     if(!fPHOSPhotonWeightFormula) 
130       fPHOSPhotonWeightFormula = new TFormula("phos_photon_weight",
131                                               fPHOSPhotonWeightFormulaExpression);
132     return fPHOSPhotonWeightFormula                                           ; } 
133   
134   TFormula * GetPHOSPi0WeightFormula()         { 
135     if(!fPHOSPi0WeightFormula) 
136       fPHOSPi0WeightFormula = new TFormula("phos_pi0_weight",
137                                            fPHOSPi0WeightFormulaExpression);
138     return fPHOSPi0WeightFormula                                              ; } 
139   
140   TString  GetPHOSPhotonWeightFormulaExpression() const { return fPHOSPhotonWeightFormulaExpression ; } 
141   TString  GetPHOSPi0WeightFormulaExpression()    const { return fPHOSPi0WeightFormulaExpression    ; } 
142   
143   //Weight setters
144   void    SetEMCALPhotonWeight  (Float_t  w)   { fEMCALPhotonWeight   = w     ; }
145   void    SetEMCALPi0Weight     (Float_t  w)   { fEMCALPi0Weight      = w     ; }
146   void    SetEMCALElectronWeight(Float_t  w)   { fEMCALElectronWeight = w     ; }
147   void    SetEMCALChargeWeight  (Float_t  w)   { fEMCALChargeWeight   = w     ; }
148   void    SetEMCALNeutralWeight (Float_t  w)   { fEMCALNeutralWeight  = w     ; }
149   void    SetPHOSPhotonWeight   (Float_t  w)   { fPHOSPhotonWeight    = w     ; }
150   void    SetPHOSPi0Weight      (Float_t  w)   { fPHOSPi0Weight       = w     ; }
151   void    SetPHOSElectronWeight (Float_t  w)   { fPHOSElectronWeight  = w     ; }
152   void    SetPHOSChargeWeight   (Float_t  w)   { fPHOSChargeWeight    = w     ; }
153   void    SetPHOSNeutralWeight  (Float_t  w)   { fPHOSNeutralWeight   = w     ; }
154   
155   void    UsePHOSPIDWeightFormula   (Bool_t ok )           { fPHOSWeightFormula                 = ok ; } 
156   void    SetPHOSPhotonWeightFormulaExpression(TString ph) { fPHOSPhotonWeightFormulaExpression = ph ; } 
157   void    SetPHOSPi0WeightFormulaExpression   (TString pi) { fPHOSPi0WeightFormulaExpression    = pi ; }
158   
159   //PID cuts 
160   
161   void    SetEMCALLambda0CutMax(Float_t lcut ) { fEMCALL0CutMax     = lcut    ; }
162   Float_t GetEMCALLambda0CutMax()        const { return fEMCALL0CutMax        ; }   
163   
164   void    SetEMCALLambda0CutMin(Float_t lcut ) { fEMCALL0CutMin     = lcut    ; }
165   Float_t GetEMCALLambda0CutMin()        const { return fEMCALL0CutMin        ; }   
166   
167   void    SetEMCALDEtaCut(Float_t dcut )       { fEMCALDEtaCut      = dcut    ; }
168   Float_t GetEMCALDEtaCut()              const { return fEMCALDEtaCut         ; }   
169   
170   void    SetEMCALDPhiCut(Float_t dcut )       { fEMCALDPhiCut      = dcut    ; }
171   Float_t GetEMCALDPhiCut()              const { return fEMCALDPhiCut         ; }   
172   
173   void    SetTOFCut(Float_t tcut )             { fTOFCut            = tcut    ; }
174   Float_t GetTOFCut()                    const { return fTOFCut               ; }   
175   
176   void    SetPHOSRCut(Float_t rcut )           { fPHOSRCut          = rcut    ; }
177   Float_t GetPHOSRCut()                  const { return fPHOSRCut             ; }   
178
179   void    SetPHOSDispersionCut(Float_t dcut )  { fPHOSDispersionCut = dcut    ; }
180   Float_t GetPHOSDispersionCut()         const { return fPHOSDispersionCut    ; }   
181   
182   
183 private:
184   
185   Int_t     fDebug;                             // Debug level
186   Int_t     fParticleFlux;                      // Particle flux for setting PID parameters
187
188   // Bayesian
189   AliEMCALPIDUtils * fEMCALPIDUtils;            // Pointer to EMCALPID to redo the PID Bayesian calculation
190   Bool_t    fUseBayesianWeights;                // Select clusters based on weights calculated in reconstruction
191   Bool_t    fRecalculateBayesian;               // Recalculate PID bayesian or use simple PID?
192
193   Float_t   fEMCALPhotonWeight;                 // Bayesian PID weight for photons in EMCAL 
194   Float_t   fEMCALPi0Weight;                    // Bayesian PID weight for pi0 in EMCAL 
195   Float_t   fEMCALElectronWeight;               // Bayesian PID weight for electrons in EMCAL 
196   Float_t   fEMCALChargeWeight;                 // Bayesian PID weight for charged hadrons in EMCAL 
197   Float_t   fEMCALNeutralWeight;                // Bayesian PID weight for neutral hadrons in EMCAL 
198   Float_t   fPHOSPhotonWeight;                  // Bayesian PID weight for photons in PHOS 
199   Float_t   fPHOSPi0Weight;                     // Bayesian PID weight for pi0 in PHOS 
200   Float_t   fPHOSElectronWeight;                // Bayesian PID weight for electrons in PHOS 
201   Float_t   fPHOSChargeWeight;                  // Bayesian PID weight for charged hadrons in PHOS 
202   Float_t   fPHOSNeutralWeight;                 // Bayesian PID weight for neutral hadrons in PHOS 
203   
204   Bool_t    fPHOSWeightFormula ;                // Use parametrized weight threshold, function of energy
205   TFormula *fPHOSPhotonWeightFormula ;          // Formula for photon weight
206   TFormula *fPHOSPi0WeightFormula ;             // Formula for pi0 weight
207   TString   fPHOSPhotonWeightFormulaExpression; // Photon weight formula in string
208   TString   fPHOSPi0WeightFormulaExpression;    // Pi0 weight formula in string
209
210   // PID calculation
211   Float_t   fEMCALL0CutMax;                     // Max Cut on shower shape lambda0, used in PID evaluation, only EMCAL
212   Float_t   fEMCALL0CutMin;                     // Min Cut on shower shape lambda0, used in PID evaluation, only EMCAL
213   Float_t   fEMCALDEtaCut;                      // Track matching cut on Dz
214   Float_t   fEMCALDPhiCut;                      // Track matching cut on Dx
215
216   Float_t   fTOFCut;                            // Cut on TOF, used in PID evaluation
217   
218   Float_t   fPHOSDispersionCut;                 // Shower shape elipse radious cut
219   Float_t   fPHOSRCut;                          // Track-Cluster distance cut for track matching in PHOS  
220   
221   AliCaloPID & operator = (const AliCaloPID & g) ; // cpy assignment
222   AliCaloPID(              const AliCaloPID & g) ; // cpy ctor
223   
224   ClassDef(AliCaloPID,11)
225 } ;
226
227
228 #endif //ALICALOPID_H
229
230
231