CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ElectronMVAEstimatorRun2Spring15Trig.h
Go to the documentation of this file.
1 #ifndef RecoEgamma_ElectronIdentification_ElectronMVAEstimatorRun2Spring15Trig_H
2 #define RecoEgamma_ElectronIdentification_ElectronMVAEstimatorRun2Spring15Trig_H
3 
5 
7 
9 
13 
15 
16 #include <vector>
17 #include <string>
18 #include <TROOT.h>
19 #include "TMVA/Factory.h"
20 #include "TMVA/Tools.h"
21 #include "TMVA/Reader.h"
22 
24 
25  public:
26 
27  // Define here the number and the meaning of the categories
28  // for this specific MVA
29  const int nCategories = 3;
31  UNDEFINED = -1,
32  CAT_EB1 = 0,
33  CAT_EB2 = 1,
34  CAT_EE = 2
35  };
36 
37  // Define the struct that contains all necessary for MVA variables
38  // Note: all variables have to be floats for TMVA Reader, even if
39  // the training was done with ints.
40  struct AllVariables {
41  // Pure ECAL -> shower shapes
42  float see; // 0
43  float spp; // 1
44  float OneMinusE1x5E5x5; // 2
45  float R9; // 3
46  float etawidth; // 4
47  float phiwidth; // 5
48  float HoE; // 6
49 
50  // Endcap only variables
51  float PreShowerOverRaw; // 7
52 
53  //Pure tracking variables
54  float kfhits; // 8
55  float kfchi2; // 9
56  float gsfchi2; // 10
57  // Energy matching
58  float fbrem; // 11
59 
60  float gsfhits; // 12
62  float convVtxFitProbability; // 14
63 
64  float EoP; // 15
65  float eleEoPout; // 16
66  float IoEmIoP; // 17
67  // Geometrical matchings
68  float deta; // 18
69  float dphi; // 19
70  float detacalo; // 20
71 
72  // Spectator variables
73  // ... none in this version ...
74 
75  // Other variables. These ones are not needed for this version
76  // of MVA, but kept in this data structure for convenience.
77  float pt; // 21
78  float SCeta; //22
79 
80 
81  };
82 
83  // Constructor and destructor
86 
87  // Calculation of the MVA value
88  float mvaValue( const edm::Ptr<reco::Candidate>& particle, const edm::Event&) const override;
89 
90  // Utility functions
91  std::unique_ptr<const GBRForest> createSingleReader(const int iCategory,
92  const edm::FileInPath &weightFile);
93 
94  virtual int getNCategories() const override { return nCategories; }
95  bool isEndcapCategory( int category ) const;
96  virtual const std::string& getName() const override final { return _name; }
97  virtual const std::string& getTag() const override final { return _tag; }
98 
99  // Functions that should work on both pat and reco electrons
100  // (use the fact that pat::Electron inherits from reco::GsfElectron)
101  std::vector<float> fillMVAVariables(const edm::Ptr<reco::Candidate>& particle, const edm::Event&) const override;
102  int findCategory( const edm::Ptr<reco::Candidate>& particle) const override;
103  // The function below ensures that the variables passed to MVA are
104  // within reasonable bounds
105  void constrainMVAVariables(AllVariables&) const;
106 
107  // Call this function once after the constructor to declare
108  // the needed event content pieces to the framework
109  void setConsumes(edm::ConsumesCollector&&) const override final;
110  // Call this function once per event to retrieve all needed
111  // event content pices
112 
113  private:
114 
115  // MVA name. This is a unique name for this MVA implementation.
116  // It will be used as part of ValueMap names.
117  // For simplicity, keep it set to the class name.
119  // MVA tag. This is an additional string variable to distinguish
120  // instances of the estimator of this class configured with different
121  // weight files.
122  const std::string _tag;
123 
124  // Data members
125  std::vector< std::unique_ptr<const GBRForest> > _gbrForests;
126 
127  // All variables needed by this MVA
128  const std::string _MethodName;
130 
131  //
132  // Declare all tokens that will be needed to retrieve misc
133  // data from the event content required by this MVA
134  //
136  // Conversions in AOD and miniAOD have different names
139 
140 
141 };
142 
143 #endif
std::unique_ptr< const GBRForest > createSingleReader(const int iCategory, const edm::FileInPath &weightFile)
std::vector< std::unique_ptr< const GBRForest > > _gbrForests
std::vector< float > fillMVAVariables(const edm::Ptr< reco::Candidate > &particle, const edm::Event &) const override
void setConsumes(edm::ConsumesCollector &&) const overridefinal
float mvaValue(const edm::Ptr< reco::Candidate > &particle, const edm::Event &) const override
virtual const std::string & getTag() const overridefinal
virtual const std::string & getName() const overridefinal
int findCategory(const edm::Ptr< reco::Candidate > &particle) const override
ElectronMVAEstimatorRun2Spring15Trig(const edm::ParameterSet &conf)
string const
Definition: compareJSON.py:14