CMS 3D CMS Logo

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