CMS 3D CMS Logo

ElectronMVAEstimator.cc
Go to the documentation of this file.
6 
8  cfg_{}
9 {}
10 
12  cfg_{}
13 {
14  // Taken from Daniele (his mail from the 30/11)
15  // tmvaReader.BookMVA("BDTSimpleCat","../Training/weights_Root527b_3Depth_DanVarConvRej_2PtBins_10Pt_800TPrune5_Min100Events_NoBjets_half/TMVA_BDTSimpleCat.weights.xm");
16  // training of the 7/12 with Nvtx added
18 }
19 
21 {
22  for(auto& weightsfile : cfg_.vweightsfiles) {
23  gbr_.push_back( GBRForestTools::createGBRForest( weightsfile ));
24  }
25 }
26 
27 double ElectronMVAEstimator::mva(const reco::GsfElectron& myElectron, int nvertices ) const {
28  float vars[18];
29 
30  vars[0] = myElectron.fbrem();
31  vars[1] = std::abs(myElectron.deltaEtaSuperClusterTrackAtVtx());
32  vars[2] = std::abs(myElectron.deltaPhiSuperClusterTrackAtVtx());
33  vars[3] = myElectron.sigmaIetaIeta();
34  vars[4] = myElectron.hcalOverEcal();
35  vars[5] = myElectron.eSuperClusterOverP();
36  vars[6] = (myElectron.e5x5()) !=0. ? 1.-(myElectron.e1x5()/myElectron.e5x5()) : -1. ;
37  vars[7] = myElectron.eEleClusterOverPout();
38  vars[8] = std::abs(myElectron.deltaEtaEleClusterTrackAtCalo());
39 
40  bool validKF= false;
41 
42  reco::TrackRef myTrackRef = myElectron.closestCtfTrackRef();
43  validKF = (myTrackRef.isNonnull() && myTrackRef.isAvailable());
44 
45  vars[9] = (validKF) ? myTrackRef->normalizedChi2() : 0 ;
46  vars[10] = (validKF) ? myTrackRef->hitPattern().trackerLayersWithMeasurement() : -1.;
47  vars[11] = myElectron.gsfTrack()->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS);
48  vars[12] = std::abs(myElectron.convDist());
49  vars[13] = std::abs(myElectron.convDcot());
50  vars[14] = nvertices;
51  vars[15] = myElectron.eta();
52  vars[16] = myElectron.pt();
53  vars[17] = myElectron.ecalDrivenSeed();
54 
55  bindVariables(vars);
56 
57  //0 pt < 10 && abs(eta)<=1.485
58  //1 pt >= 10 && abs(eta)<=1.485
59  //2 pt < 10 && abs(eta)> 1.485
60  //3 pt >= 10 && abs(eta)> 1.485
61 
62  const unsigned index = (unsigned)(myElectron.pt() >= 10) + 2*(unsigned)(std::abs(myElectron.eta()) > 1.485);
63 
64  double result = gbr_[index]->GetAdaBoostClassifier(vars);
65 //
66 // std::cout << "fbrem" << vars[0] << std::endl;
67 // std::cout << "detain"<< vars[1] << std::endl;
68 // std::cout << "dphiin"<< vars[2] << std::endl;
69 // std::cout << "sieie"<< vars[3] << std::endl;
70 // std::cout << "hoe"<< vars[4] << std::endl;
71 // std::cout << "eop"<< vars[5] << std::endl;
72 // std::cout << "e1x5e5x5"<< vars[6] << std::endl;
73 // std::cout << "eleopout"<< vars[7] << std::endl;
74 // std::cout << "detaeleout"<< vars[8] << std::endl;
75 // std::cout << "kfchi2"<< vars[9] << std::endl;
76 // std::cout << "kfhits"<< vars[10] << std::endl;
77 // std::cout << "mishits"<<vars[11] << std::endl;
78 // std::cout << "dist"<< vars[12] << std::endl;
79 // std::cout << "dcot"<< vars[13] << std::endl;
80 // std::cout << "nvtx"<< vars[14] << std::endl;
81 // std::cout << "eta"<< vars[15] << std::endl;
82 // std::cout << "pt"<< vars[16] << std::endl;
83 // std::cout << "ecalseed"<< vars[17] << std::endl;
84 //
85 // std::cout << " MVA " << result << std::endl;
86  return result;
87 }
88 
89 
90 void ElectronMVAEstimator::bindVariables(float vars[18]) const {
91  if(vars[0] < -1.)
92  vars[1] = -1.;
93 
94  if(vars[1] > 0.06)
95  vars[1] = 0.06;
96 
97  if(vars[2] > 0.6)
98  vars[2] = 0.6;
99 
100  if(vars[5] > 20.)
101  vars[5] = 20.;
102 
103  if(vars[7] > 20.)
104  vars[7] = 20;
105 
106  if(vars[8] > 0.2)
107  vars[8] = 0.2;
108 
109  if(vars[9] < 0.)
110  vars[9] = 0.;
111 
112  if(vars[9] > 15.)
113  vars[9] = 15.;
114 
115  if(vars[6] < -1.)
116  vars[6] = -1;
117 
118  if(vars[6] > 2.)
119  vars[6] = 2.;
120 
121  if(vars[12] > 15.)
122  vars[12] = 15.;
123 
124  if(vars[13] > 3.)
125  vars[13] = 3.;
126 }
bool isAvailable() const
Definition: Ref.h:577
GsfTrackRef gsfTrack() const override
reference to a GsfTrack
Definition: GsfElectron.h:185
double mva(const reco::GsfElectron &myElectron, int nvertices=0) const
virtual TrackRef closestCtfTrackRef() const
Definition: GsfElectron.h:201
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:253
double eta() const final
momentum pseudorapidity
float eSuperClusterOverP() const
Definition: GsfElectron.h:245
std::vector< std::unique_ptr< const GBRForest > > gbr_
double pt() const final
transverse momentum
void bindVariables(float vars[18]) const
float fbrem() const
Definition: GsfElectron.h:761
float convDist() const
Definition: GsfElectron.h:609
float convDcot() const
Definition: GsfElectron.h:610
float deltaEtaSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:249
float sigmaIetaIeta() const
Definition: GsfElectron.h:436
static std::unique_ptr< const GBRForest > createGBRForest(const std::string &weightFile)
float deltaPhiSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:252
const Configuration cfg_
float hcalOverEcal() const
Definition: GsfElectron.h:444
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float eEleClusterOverPout() const
Definition: GsfElectron.h:248
float deltaEtaEleClusterTrackAtCalo() const
Definition: GsfElectron.h:251
float e1x5() const
Definition: GsfElectron.h:438
std::vector< std::string > vweightsfiles
float e5x5() const
Definition: GsfElectron.h:440
bool ecalDrivenSeed() const
Definition: GsfElectron.h:188