test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ElectronMVAEstimator.cc
Go to the documentation of this file.
7 
8 #include "TMVA/Reader.h"
9 #include "TMVA/MethodBDT.h"
10 #include "TMVA/MethodCategory.h"
11 
12 namespace {
13  constexpr char ele_mva_name[] = "BDTSimpleCat";
14 }
15 
17  cfg_{}
18 {}
19 
21  cfg_{}
22 {
23  TMVA::Reader tmvaReader("!Color:Silent");
24  tmvaReader.AddVariable("fbrem",&fbrem);
25  tmvaReader.AddVariable("detain", &detain);
26  tmvaReader.AddVariable("dphiin", &dphiin);
27  tmvaReader.AddVariable("sieie", &sieie);
28  tmvaReader.AddVariable("hoe", &hoe);
29  tmvaReader.AddVariable("eop", &eop);
30  tmvaReader.AddVariable("e1x5e5x5", &e1x5e5x5);
31  tmvaReader.AddVariable("eleopout", &eleopout);
32  tmvaReader.AddVariable("detaeleout", &detaeleout);
33  tmvaReader.AddVariable("kfchi2", &kfchi2);
34  tmvaReader.AddVariable("kfhits", &mykfhits);
35  tmvaReader.AddVariable("mishits",&mymishits);
36  tmvaReader.AddVariable("dist", &absdist);
37  tmvaReader.AddVariable("dcot", &absdcot);
38  tmvaReader.AddVariable("nvtx", &myNvtx);
39 
40  tmvaReader.AddSpectator("eta",&eta);
41  tmvaReader.AddSpectator("pt",&pt);
42  tmvaReader.AddSpectator("ecalseed",&ecalseed);
43 
44  // Taken from Daniele (his mail from the 30/11)
45  // tmvaReader.BookMVA("BDTSimpleCat","../Training/weights_Root527b_3Depth_DanVarConvRej_2PtBins_10Pt_800TPrune5_Min100Events_NoBjets_half/TMVA_BDTSimpleCat.weights.xm");
46  // training of the 7/12 with Nvtx added
47  std::unique_ptr<TMVA::IMethod> temp( tmvaReader.BookMVA(ele_mva_name,fileName.c_str()) );
48  gbr.emplace_back(new GBRForest( dynamic_cast<TMVA::MethodBDT*>( tmvaReader.FindMVA(ele_mva_name) ) ) );
49 }
50 
52  std::vector<std::string> weightsfiles;
53  std::string path_mvaWeightFileEleID;
54  for(unsigned ifile=0 ; ifile < cfg_.vweightsfiles.size() ; ++ifile) {
55  path_mvaWeightFileEleID = edm::FileInPath ( cfg_.vweightsfiles[ifile].c_str() ).fullPath();
56  weightsfiles.push_back(path_mvaWeightFileEleID);
57  }
58  for( const auto& wgtfile : weightsfiles ) {
59  TMVA::Reader tmvaReader("!Color:Silent");
60  tmvaReader.AddVariable("fbrem",&fbrem);
61  tmvaReader.AddVariable("detain", &detain);
62  tmvaReader.AddVariable("dphiin", &dphiin);
63  tmvaReader.AddVariable("sieie", &sieie);
64  tmvaReader.AddVariable("hoe", &hoe);
65  tmvaReader.AddVariable("eop", &eop);
66  tmvaReader.AddVariable("e1x5e5x5", &e1x5e5x5);
67  tmvaReader.AddVariable("eleopout", &eleopout);
68  tmvaReader.AddVariable("detaeleout", &detaeleout);
69  tmvaReader.AddVariable("kfchi2", &kfchi2);
70  tmvaReader.AddVariable("kfhits", &mykfhits);
71  tmvaReader.AddVariable("mishits",&mymishits);
72  tmvaReader.AddVariable("dist", &absdist);
73  tmvaReader.AddVariable("dcot", &absdcot);
74  tmvaReader.AddVariable("nvtx", &myNvtx);
75 
76  tmvaReader.AddSpectator("eta",&eta);
77  tmvaReader.AddSpectator("pt",&pt);
78  tmvaReader.AddSpectator("ecalseed",&ecalseed);
79 
80  // Taken from Daniele (his mail from the 30/11)
81  // tmvaReader.BookMVA("BDTSimpleCat","../Training/weights_Root527b_3Depth_DanVarConvRej_2PtBins_10Pt_800TPrune5_Min100Events_NoBjets_half/TMVA_BDTSimpleCat.weights.xm");
82  // training of the 7/12 with Nvtx added
83  std::unique_ptr<TMVA::IMethod> temp( tmvaReader.BookMVA(ele_mva_name,wgtfile) );
84  gbr.emplace_back(new GBRForest( dynamic_cast<TMVA::MethodBDT*>( tmvaReader.FindMVA(ele_mva_name) ) ) );
85  }
86 }
87 
88 double ElectronMVAEstimator::mva(const reco::GsfElectron& myElectron, int nvertices ) const {
89  float vars[18];
90 
91  vars[0] = myElectron.fbrem();
92  vars[1] = std::abs(myElectron.deltaEtaSuperClusterTrackAtVtx());
93  vars[2] = std::abs(myElectron.deltaPhiSuperClusterTrackAtVtx());
94  vars[3] = myElectron.sigmaIetaIeta();
95  vars[4] = myElectron.hcalOverEcal();
96  vars[5] = myElectron.eSuperClusterOverP();
97  vars[6] = (myElectron.e5x5()) !=0. ? 1.-(myElectron.e1x5()/myElectron.e5x5()) : -1. ;
98  vars[7] = myElectron.eEleClusterOverPout();
99  vars[8] = std::abs(myElectron.deltaEtaEleClusterTrackAtCalo());
100 
101  bool validKF= false;
102 
103  reco::TrackRef myTrackRef = myElectron.closestCtfTrackRef();
104  validKF = (myTrackRef.isNonnull() && myTrackRef.isAvailable());
105 
106  vars[9] = (validKF) ? myTrackRef->normalizedChi2() : 0 ;
107  vars[10] = (validKF) ? myTrackRef->hitPattern().trackerLayersWithMeasurement() : -1.;
108  vars[11] = myElectron.gsfTrack()->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS);
109  vars[12] = std::abs(myElectron.convDist());
110  vars[13] = std::abs(myElectron.convDcot());
111  vars[14] = nvertices;
112  vars[15] = myElectron.eta();
113  vars[16] = myElectron.pt();
114  vars[17] = myElectron.ecalDrivenSeed();
115 
116  bindVariables(vars);
117 
118  //0 pt &lt; 10 &amp;&amp; abs(eta)&lt;=1.485
119  //1 pt &gt;= 10 &amp;&amp; abs(eta)&lt;=1.485
120  //2 pt &lt; 10 &amp;&amp; abs(eta)&gt; 1.485
121  //3 pt &gt;= 10 &amp;&amp; abs(eta)&gt; 1.485
122 
123  const unsigned index = (unsigned)(myElectron.pt() >= 10) + 2*(unsigned)(std::abs(myElectron.eta()) > 1.485);
124 
125  double result = gbr[index]->GetAdaBoostClassifier(vars);
126 //
127 // std::cout << "fbrem" << vars[0] << std::endl;
128 // std::cout << "detain"<< vars[1] << std::endl;
129 // std::cout << "dphiin"<< vars[2] << std::endl;
130 // std::cout << "sieie"<< vars[3] << std::endl;
131 // std::cout << "hoe"<< vars[4] << std::endl;
132 // std::cout << "eop"<< vars[5] << std::endl;
133 // std::cout << "e1x5e5x5"<< vars[6] << std::endl;
134 // std::cout << "eleopout"<< vars[7] << std::endl;
135 // std::cout << "detaeleout"<< vars[8] << std::endl;
136 // std::cout << "kfchi2"<< vars[9] << std::endl;
137 // std::cout << "kfhits"<< vars[10] << std::endl;
138 // std::cout << "mishits"<<vars[11] << std::endl;
139 // std::cout << "dist"<< vars[12] << std::endl;
140 // std::cout << "dcot"<< vars[13] << std::endl;
141 // std::cout << "nvtx"<< vars[14] << std::endl;
142 // std::cout << "eta"<< vars[15] << std::endl;
143 // std::cout << "pt"<< vars[16] << std::endl;
144 // std::cout << "ecalseed"<< vars[17] << std::endl;
145 //
146 // std::cout << " MVA " << result << std::endl;
147  return result;
148 }
149 
150 
151 void ElectronMVAEstimator::bindVariables(float vars[18]) const {
152  if(vars[0] < -1.)
153  vars[1] = -1.;
154 
155  if(vars[1] > 0.06)
156  vars[1] = 0.06;
157 
158  if(vars[2] > 0.6)
159  vars[2] = 0.6;
160 
161  if(vars[5] > 20.)
162  vars[5] = 20.;
163 
164  if(vars[7] > 20.)
165  vars[7] = 20;
166 
167  if(vars[8] > 0.2)
168  vars[8] = 0.2;
169 
170  if(vars[9] < 0.)
171  vars[9] = 0.;
172 
173  if(vars[9] > 15.)
174  vars[9] = 15.;
175 
176  if(vars[6] < -1.)
177  vars[6] = -1;
178 
179  if(vars[6] > 2.)
180  vars[6] = 2.;
181 
182  if(vars[12] > 15.)
183  vars[12] = 15.;
184 
185  if(vars[13] > 3.)
186  vars[13] = 3.;
187 }
bool isAvailable() const
Definition: Ref.h:576
double mva(const reco::GsfElectron &myElectron, int nvertices=0) const
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:252
tuple cfg
Definition: looper.py:293
float eSuperClusterOverP() const
Definition: GsfElectron.h:243
void bindVariables(float vars[18]) const
float fbrem() const
Definition: GsfElectron.h:698
float convDist() const
Definition: GsfElectron.h:553
#define constexpr
TrackRef closestCtfTrackRef() const
Definition: GsfElectron.h:199
tuple result
Definition: mps_fire.py:84
float convDcot() const
Definition: GsfElectron.h:554
float deltaEtaSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:247
float sigmaIetaIeta() const
Definition: GsfElectron.h:416
float deltaPhiSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:250
const Configuration cfg_
float hcalOverEcal() const
Definition: GsfElectron.h:424
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float eEleClusterOverPout() const
Definition: GsfElectron.h:246
float deltaEtaEleClusterTrackAtCalo() const
Definition: GsfElectron.h:249
float e1x5() const
Definition: GsfElectron.h:418
std::vector< std::string > vweightsfiles
float e5x5() const
Definition: GsfElectron.h:420
virtual double eta() const final
momentum pseudorapidity
virtual GsfTrackRef gsfTrack() const
reference to a GsfTrack
Definition: GsfElectron.h:183
virtual double pt() const final
transverse momentum
bool ecalDrivenSeed() const
Definition: GsfElectron.h:186
std::vector< std::unique_ptr< const GBRForest > > gbr