]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/CaloTrackCorrBase/AliCaloPID.h
Move cluster splitting method from AliAnaInsideClusterInvariantMass to AliCalorimeter...
[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 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 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 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 AliVCaloCells;
41 class AliAODPWG4Particle;
42 class AliEMCALPIDUtils;
43 class AliCalorimeterUtils;
44 class AliVEvent;
45
46 class AliCaloPID : public TObject {
47         
48  public: 
49   
50   AliCaloPID() ; // ctor
51   AliCaloPID(const Int_t particleFlux) ; // ctor, to be used when recalculating bayesian PID
52   AliCaloPID(const TNamed * emcalpid) ; // ctor, to be used when recalculating bayesian PID and need different parameters
53   virtual ~AliCaloPID() ;//virtual dtor
54         
55   enum PidType 
56   {
57     kPhoton         = 22,
58     kPi0            = 111,
59     kEta            = 221, 
60     kElectron       = 11, 
61     kEleCon         =-11, 
62     kNeutralHadron  = 2112, 
63     kChargedHadron  = 211, 
64     kNeutralUnknown = 130, 
65     kChargedUnknown = 321
66   };
67   
68   enum TagType {kPi0Decay, kEtaDecay, kOtherDecay, kConversion, kNoTag = -1};
69   
70   // Main methods
71   
72   TList *   GetCreateOutputObjects();
73
74   void      InitParameters();
75   
76   Int_t     GetIdentifiedParticleTypeFromBayesWeights(const Bool_t isEMCAL, const Double_t * pid, const Float_t energy) ;
77
78   Int_t     GetIdentifiedParticleTypeFromClusterSplitting(AliVCluster * cluster, AliVCaloCells* cells, 
79                                                           AliCalorimeterUtils * caloutils,
80                                                           Double_t vertex[3], 
81                                                           Int_t & nLocMax, Double_t & mass, Double_t & angle) ;
82   
83   Int_t     GetIdentifiedParticleType(const AliVCluster * cluster) ;
84   
85   TString   GetPIDParametersList();
86   
87   Bool_t    IsTrackMatched(AliVCluster * cluster, AliCalorimeterUtils* cu, AliVEvent* event) const ;    
88   
89   void      SetPIDBits(AliVCluster * cluster, AliAODPWG4Particle *aodph, 
90                        AliCalorimeterUtils* cu, AliVEvent* event);
91   
92   void      Print(const Option_t * opt)const;
93   
94   void      PrintClusterPIDWeights(const Double_t * pid) const;
95   
96   //Check if cluster photon-like. Uses photon cluster parameterization in real pp data 
97   //Returns distance in sigmas. Recommended cut 2.5
98   Float_t TestPHOSDispersion(const Double_t pt, const Double_t m20, const Double_t m02) const ; 
99   //Checks distance to the closest track. Takes into account 
100   //non-perpendicular incidence of tracks.
101   Float_t TestPHOSChargedVeto(const Double_t dx,  const Double_t dz, const Double_t ptTrack, 
102                               const Int_t chargeTrack, const Double_t mf) const ;
103   
104   // Setters, getters
105   
106   void    SetDebug(Int_t deb)                  { fDebug = deb                 ; }
107   Int_t   GetDebug()                     const { return fDebug                ; }       
108
109   enum    eventType{kLow,kHigh};
110   void    SetLowParticleFlux()                 { fParticleFlux        = kLow  ; }
111   void    SetHighParticleFlux()                { fParticleFlux        = kHigh ; }  
112   // not really used, only for bayesian recalculation in EMCAL, but could be useful in future
113   
114   // Bayesian
115   
116   void    SwitchOnBayesian()                   { fUseBayesianWeights  = kTRUE ; }
117   void    SwitchOffBayesian()                  { fUseBayesianWeights  = kFALSE; }
118   void    SwitchOnBayesianRecalculation()      { fRecalculateBayesian = kTRUE ; fUseBayesianWeights  = kTRUE ;} // EMCAL
119   void    SwitchOffBayesianRecalculation()     { fRecalculateBayesian = kFALSE; } // EMCAL
120   
121   AliEMCALPIDUtils * GetEMCALPIDUtils() ; 
122   
123   //Weight getters
124   Float_t GetEMCALPhotonWeight()         const { return fEMCALPhotonWeight    ; }
125   Float_t GetEMCALPi0Weight()            const { return fEMCALPi0Weight       ; }
126   Float_t GetEMCALElectronWeight()       const { return fEMCALElectronWeight  ; }
127   Float_t GetEMCALChargeWeight()         const { return fEMCALChargeWeight    ; }
128   Float_t GetEMCALNeutralWeight()        const { return fEMCALNeutralWeight   ; }
129   Float_t GetPHOSPhotonWeight()          const { return fPHOSPhotonWeight     ; }
130   Float_t GetPHOSPi0Weight()             const { return fPHOSPi0Weight        ; }
131   Float_t GetPHOSElectronWeight()        const { return fPHOSElectronWeight   ; }
132   Float_t GetPHOSChargeWeight()          const { return fPHOSChargeWeight     ; }
133   Float_t GetPHOSNeutralWeight()         const { return fPHOSNeutralWeight    ; }
134   
135   Bool_t  IsPHOSPIDWeightFormulaOn()     const { return fPHOSWeightFormula    ; } 
136
137   TFormula * GetPHOSPhotonWeightFormula()      { 
138     if(!fPHOSPhotonWeightFormula) 
139       fPHOSPhotonWeightFormula = new TFormula("phos_photon_weight",
140                                               fPHOSPhotonWeightFormulaExpression);
141     return fPHOSPhotonWeightFormula                                           ; } 
142   
143   TFormula * GetPHOSPi0WeightFormula()         { 
144     if(!fPHOSPi0WeightFormula) 
145       fPHOSPi0WeightFormula = new TFormula("phos_pi0_weight",
146                                            fPHOSPi0WeightFormulaExpression);
147     return fPHOSPi0WeightFormula                                              ; } 
148   
149   TString  GetPHOSPhotonWeightFormulaExpression() const { return fPHOSPhotonWeightFormulaExpression ; } 
150   TString  GetPHOSPi0WeightFormulaExpression()    const { return fPHOSPi0WeightFormulaExpression    ; } 
151   
152   //Weight setters
153   void    SetEMCALPhotonWeight  (Float_t  w)   { fEMCALPhotonWeight   = w     ; }
154   void    SetEMCALPi0Weight     (Float_t  w)   { fEMCALPi0Weight      = w     ; }
155   void    SetEMCALElectronWeight(Float_t  w)   { fEMCALElectronWeight = w     ; }
156   void    SetEMCALChargeWeight  (Float_t  w)   { fEMCALChargeWeight   = w     ; }
157   void    SetEMCALNeutralWeight (Float_t  w)   { fEMCALNeutralWeight  = w     ; }
158   void    SetPHOSPhotonWeight   (Float_t  w)   { fPHOSPhotonWeight    = w     ; }
159   void    SetPHOSPi0Weight      (Float_t  w)   { fPHOSPi0Weight       = w     ; }
160   void    SetPHOSElectronWeight (Float_t  w)   { fPHOSElectronWeight  = w     ; }
161   void    SetPHOSChargeWeight   (Float_t  w)   { fPHOSChargeWeight    = w     ; }
162   void    SetPHOSNeutralWeight  (Float_t  w)   { fPHOSNeutralWeight   = w     ; }
163   
164   void    UsePHOSPIDWeightFormula   (Bool_t ok )           { fPHOSWeightFormula                 = ok ; } 
165   void    SetPHOSPhotonWeightFormulaExpression(TString ph) { fPHOSPhotonWeightFormulaExpression = ph ; } 
166   void    SetPHOSPi0WeightFormulaExpression   (TString pi) { fPHOSPi0WeightFormulaExpression    = pi ; }
167   
168   //PID cuts 
169   
170   void    SetEMCALLambda0CutMax(Float_t lcut ) { fEMCALL0CutMax     = lcut    ; }
171   Float_t GetEMCALLambda0CutMax()        const { return fEMCALL0CutMax        ; }   
172   
173   void    SetEMCALLambda0CutMin(Float_t lcut ) { fEMCALL0CutMin     = lcut    ; }
174   Float_t GetEMCALLambda0CutMin()        const { return fEMCALL0CutMin        ; }   
175   
176   void    SetEMCALDEtaCut(Float_t dcut )       { fEMCALDEtaCut      = dcut    ; }
177   Float_t GetEMCALDEtaCut()              const { return fEMCALDEtaCut         ; }   
178   
179   void    SetEMCALDPhiCut(Float_t dcut )       { fEMCALDPhiCut      = dcut    ; }
180   Float_t GetEMCALDPhiCut()              const { return fEMCALDPhiCut         ; }   
181   
182   void    SetTOFCut(Float_t tcut )             { fTOFCut            = tcut    ; }
183   Float_t GetTOFCut()                    const { return fTOFCut               ; }   
184   
185   void    SetPHOSRCut(Float_t rcut )           { fPHOSRCut          = rcut    ; }
186   Float_t GetPHOSRCut()                  const { return fPHOSRCut             ; }   
187
188   void    SetPHOSDispersionCut(Float_t dcut )  { fPHOSDispersionCut = dcut    ; }
189   Float_t GetPHOSDispersionCut()         const { return fPHOSDispersionCut    ; }   
190   
191   // Cluster splitting analysis
192   
193   void    SwitchOnClusterSplittingPID()        { fDoClusterSplitting = kTRUE  ; }
194   void    SwitchOffClusterplittingPID()        { fDoClusterSplitting = kFALSE ; }
195
196   void    SetClusterSplittingM02Cut(Float_t min=0, Float_t max=100) 
197   { fSplitM02MinCut   = min ; fSplitM02MaxCut  = max ; }
198   
199   void    SetClusterSplittingMinNCells(Int_t cut)   { fSplitMinNCells   = cut ; }  
200   
201   Float_t GetPi0MinMass()                const { return fMassPi0Min           ; }
202   Float_t GetEtaMinMass()                const { return fMassEtaMin           ; }
203   Float_t GetPhotonMinMass()             const { return fMassPhoMin           ; }  
204   Float_t GetPi0MaxMass()                const { return fMassPi0Max           ; }
205   Float_t GetEtaMaxMass()                const { return fMassEtaMax           ; }
206   Float_t GetPhotonMaxMass()             const { return fMassPhoMax           ; }
207   
208   void    SetPi0MassRange(Float_t min, Float_t max)    { fMassPi0Min  = min ; fMassPi0Max = max ; }
209   void    SetEtaMassRange(Float_t min, Float_t max)    { fMassEtaMin  = min ; fMassEtaMax = max ; }
210   void    SetPhotonMassRange(Float_t min, Float_t max) { fMassPhoMin  = min ; fMassPhoMax = max ; }
211   
212 private:
213   
214   Int_t     fDebug;                             // Debug level
215   Int_t     fParticleFlux;                      // Particle flux for setting PID parameters
216
217   // Bayesian
218   AliEMCALPIDUtils * fEMCALPIDUtils;            // Pointer to EMCALPID to redo the PID Bayesian calculation
219   Bool_t    fUseBayesianWeights;                // Select clusters based on weights calculated in reconstruction
220   Bool_t    fRecalculateBayesian;               // Recalculate PID bayesian or use simple PID?
221
222   Float_t   fEMCALPhotonWeight;                 // Bayesian PID weight for photons in EMCAL 
223   Float_t   fEMCALPi0Weight;                    // Bayesian PID weight for pi0 in EMCAL 
224   Float_t   fEMCALElectronWeight;               // Bayesian PID weight for electrons in EMCAL 
225   Float_t   fEMCALChargeWeight;                 // Bayesian PID weight for charged hadrons in EMCAL 
226   Float_t   fEMCALNeutralWeight;                // Bayesian PID weight for neutral hadrons in EMCAL 
227   Float_t   fPHOSPhotonWeight;                  // Bayesian PID weight for photons in PHOS 
228   Float_t   fPHOSPi0Weight;                     // Bayesian PID weight for pi0 in PHOS 
229   Float_t   fPHOSElectronWeight;                // Bayesian PID weight for electrons in PHOS 
230   Float_t   fPHOSChargeWeight;                  // Bayesian PID weight for charged hadrons in PHOS 
231   Float_t   fPHOSNeutralWeight;                 // Bayesian PID weight for neutral hadrons in PHOS 
232   
233   Bool_t    fPHOSWeightFormula ;                // Use parametrized weight threshold, function of energy
234   TFormula *fPHOSPhotonWeightFormula ;          // Formula for photon weight
235   TFormula *fPHOSPi0WeightFormula ;             // Formula for pi0 weight
236   TString   fPHOSPhotonWeightFormulaExpression; // Photon weight formula in string
237   TString   fPHOSPi0WeightFormulaExpression;    // Pi0 weight formula in string
238
239   // PID calculation
240   Float_t   fEMCALL0CutMax;                     // Max Cut on shower shape lambda0, used in PID evaluation, only EMCAL
241   Float_t   fEMCALL0CutMin;                     // Min Cut on shower shape lambda0, used in PID evaluation, only EMCAL
242   Float_t   fEMCALDEtaCut;                      // Track matching cut on Dz
243   Float_t   fEMCALDPhiCut;                      // Track matching cut on Dx
244
245   Float_t   fTOFCut;                            // Cut on TOF, used in PID evaluation
246   
247   Float_t   fPHOSDispersionCut;                 // Shower shape elipse radious cut
248   Float_t   fPHOSRCut;                          // Track-Cluster distance cut for track matching in PHOS  
249   
250   // Cluster splitting mass ranges
251   Bool_t    fDoClusterSplitting;                // Cluster splitting analysis
252   Float_t   fSplitM02MaxCut ;                   // Study clusters with l0 smaller than cut
253   Float_t   fSplitM02MinCut ;                   // Study clusters with l0 larger than cut
254   Int_t     fSplitMinNCells ;                   // Study clusters with ncells larger than cut  
255   Float_t   fMassEtaMin  ;                      // Min Eta mass
256   Float_t   fMassEtaMax  ;                      // Max Eta mass  
257   Float_t   fMassPi0Min  ;                      // Min Pi0 mass
258   Float_t   fMassPi0Max  ;                      // Min Pi0 mass
259   Float_t   fMassPhoMin  ;                      // Min Photon mass
260   Float_t   fMassPhoMax  ;                      // Min Photon mass
261   
262   AliCaloPID & operator = (const AliCaloPID & g) ; // cpy assignment
263   AliCaloPID(              const AliCaloPID & g) ; // cpy ctor
264   
265   ClassDef(AliCaloPID,12)
266   
267 } ;
268
269
270 #endif //ALICALOPID_H
271
272
273