CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ElectronMVAEstimatorRun2Phys14NonTrig.h
Go to the documentation of this file.
1 #ifndef RecoEgamma_ElectronIdentification_ElectronMVAEstimatorRun2Phys14NonTrig_H
2 #define RecoEgamma_ElectronIdentification_ElectronMVAEstimatorRun2Phys14NonTrig_H
3 
5 
7 
9 
10 #include <vector>
11 #include <string>
12 #include <memory>
13 #include <TROOT.h>
14 #include "TMVA/Factory.h"
15 #include "TMVA/Tools.h"
16 #include "TMVA/Reader.h"
17 #include "TMVA/MethodBDT.h"
18 
20 
21  public:
22 
23  // Define here the number and the meaning of the categories
24  // for this specific MVA
25  static constexpr int nCategories = 6;
27  UNDEFINED = -1,
34  };
35 
36  // Define the struct that contains all necessary for MVA variables
37  struct AllVariables {
38  float kfhits; // 0
39  // Pure ECAL -> shower shapes
40  float see; // 1
41  float spp; // 2
42  float OneMinusE1x5E5x5; // 3
43  float R9; // 4
44  float etawidth; // 5
45  float phiwidth; // 6
46  float HoE; // 7
47  // Endcap only variables
48  float PreShowerOverRaw; // 8
49  //Pure tracking variables
50  float kfchi2; // 9
51  float gsfchi2; // 10
52  // Energy matching
53  float fbrem; // 11
54  float EoP; // 12
55  float eleEoPout; // 13
56  float IoEmIoP; // 14
57  // Geometrical matchings
58  float deta; // 15
59  float dphi; // 16
60  float detacalo; // 17
61  // Spectator variables
62  float pt; // 18
63  float isBarrel; // 19
64  float isEndcap; // 20
65  float SCeta; // 21
66  };
67 
68  // Constructor and destructor
71 
72  // Calculation of the MVA value
73  float mvaValue( const edm::Ptr<reco::Candidate>& particle, const edm::Event& evt) const;
74 
75  // Utility functions
76  std::unique_ptr<const GBRForest> createSingleReader(const int iCategory, const edm::FileInPath &weightFile) ;
77 
78  virtual int getNCategories() const override final { return nCategories; }
79  bool isEndcapCategory( int category ) const;
80  virtual const std::string& getName() const override final { return _name; }
81  virtual const std::string& getTag() const override final { return _tag; }
82 
83  // Functions that should work on both pat and reco electrons
84  // (use the fact that pat::Electron inherits from reco::GsfElectron)
85  std::vector<float> fillMVAVariables(const edm::Ptr<reco::Candidate>& particle, const edm::Event&) const;
86  int findCategory(const edm::Ptr<reco::Candidate>& particle) const;
87  // The function below ensures that the variables passed to MVA are
88  // within reasonable bounds
89  void constrainMVAVariables(AllVariables& vars) const;
90 
91  private:
92 
93  // MVA name. This is a unique name for this MVA implementation.
94  // It will be used as part of ValueMap names.
95  // For simplicity, keep it set to the class name.
96  const std::string _name = "ElectronMVAEstimatorRun2Phys14NonTrig";
97  // MVA tag. This is an additional string variable to distinguish
98  // instances of the estimator of this class configured with different
99  // weight files.
101 
102  // Data members
103  std::vector< std::unique_ptr<const GBRForest> > _gbrForests;
104 
105  // All variables needed by this MVA
108 
109 };
110 
111 #endif
std::vector< float > fillMVAVariables(const edm::Ptr< reco::Candidate > &particle, const edm::Event &) const
#define constexpr
int findCategory(const edm::Ptr< reco::Candidate > &particle) const
std::vector< std::unique_ptr< const GBRForest > > _gbrForests
std::unique_ptr< const GBRForest > createSingleReader(const int iCategory, const edm::FileInPath &weightFile)
ElectronMVAEstimatorRun2Phys14NonTrig(const edm::ParameterSet &conf)
virtual const std::string & getTag() const overridefinal
float mvaValue(const edm::Ptr< reco::Candidate > &particle, const edm::Event &evt) const
virtual const std::string & getName() const overridefinal
string const
Definition: compareJSON.py:14