CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Static Public Attributes | Private Member Functions | Private Attributes
SoftElectronMVAEstimator Class Reference

#include <SoftElectronMVAEstimator.h>

Classes

struct  Configuration
 

Public Member Functions

double mva (const reco::GsfElectron &myElectron, const reco::VertexCollection &) const
 
 SoftElectronMVAEstimator (const Configuration &)
 
 ~SoftElectronMVAEstimator ()
 

Static Public Attributes

static unsigned int ExpectedNBins = 1
 

Private Member Functions

void bindVariables (float vars[25]) const
 
void init ()
 

Private Attributes

const Configuration cfg_
 
std::vector< std::unique_ptr< const GBRForest > > gbr_
 

Detailed Description

Definition at line 14 of file SoftElectronMVAEstimator.h.

Constructor & Destructor Documentation

SoftElectronMVAEstimator::SoftElectronMVAEstimator ( const Configuration cfg)

Definition at line 8 of file SoftElectronMVAEstimator.cc.

References cfg_, createGBRForest(), ExpectedNBins, gbr_, and SoftElectronMVAEstimator::Configuration::vweightsfiles.

8  :cfg_(cfg)
9 {
10  //Check number of weight files given
11  if (ExpectedNBins != cfg_.vweightsfiles.size() ) {
12  edm::LogError ("Soft Electron MVA Error") <<
13  "Expected Number of bins = " << ExpectedNBins <<
14  " does not equal to weightsfiles.size() = " <<
15  cfg_.vweightsfiles.size() << std::endl;
16  }
17 
18  for(auto& weightsfile : cfg_.vweightsfiles) {
19  // Taken from Daniele (his mail from the 30/11)
20  // training of the 7/12 with Nvtx added
21  gbr_.push_back(createGBRForest( weightsfile ));
22  }
23 }
std::vector< std::unique_ptr< const GBRForest > > gbr_
std::unique_ptr< const GBRForest > createGBRForest(const std::string &weightsFile)
SoftElectronMVAEstimator::~SoftElectronMVAEstimator ( )

Definition at line 26 of file SoftElectronMVAEstimator.cc.

27 { }

Member Function Documentation

void SoftElectronMVAEstimator::bindVariables ( float  vars[25]) const
private

Definition at line 100 of file SoftElectronMVAEstimator.cc.

References funct::abs().

Referenced by mva().

100  {
101  if( vars[0] < -1.) //fbrem
102  vars[0] = -1.;
103 
104  vars[11] = std::abs(vars[11]); // deta
105  if(vars[11] > 0.06)
106  vars[11] = 0.06;
107 
108  vars[12] = std::abs(vars[12]);
109  if(vars[12] > 0.6)
110  vars[12] = 0.6;
111 
112  //if(EoP > 20.)
113  // EoP = 20.;
114 
115  if(vars[2] > 20.) //eleEoPout
116  vars[2] = 20.;
117 
118  vars[13] = std::abs(vars[13]); //detacalo
119  if(vars[13] > 0.2)
120  vars[13] = 0.2;
121 
122  if( vars[19] < -1.) //OneMinusE1x5E5x5
123  vars[19] = -1;
124 
125  if( vars[19] > 2.) //OneMinusE1x5E5x5
126  vars[19] = 2.;
127 
128  if( vars[7] > 200.) //gsfchi2
129  vars[7] = 200;
130 
131  if( vars[8] > 10.) //kfchi2
132  vars[8] = 10.;
133 
134 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
vars
Definition: DeepTauId.cc:77
void SoftElectronMVAEstimator::init ( )
private
double SoftElectronMVAEstimator::mva ( const reco::GsfElectron myElectron,
const reco::VertexCollection pvc 
) const

Definition at line 29 of file SoftElectronMVAEstimator.cc.

References bindVariables(), reco::GsfElectron::closestCtfTrackRef(), reco::GsfElectron::deltaEtaEleClusterTrackAtCalo(), reco::GsfElectron::deltaEtaSeedClusterTrackAtCalo(), reco::GsfElectron::deltaEtaSuperClusterTrackAtVtx(), reco::GsfElectron::deltaPhiSuperClusterTrackAtVtx(), reco::GsfElectron::e1x5(), reco::GsfElectron::e5x5(), reco::GsfElectron::ecalEnergy(), reco::GsfElectron::eEleClusterOverPout(), reco::GsfElectron::eSuperClusterOverP(), reco::LeafCandidate::eta(), reco::GsfElectron::fbrem(), gbr_, edm::Ref< C, T, F >::get(), reco::GsfElectron::gsfTrack(), reco::GsfElectron::hcalOverEcalBc(), edm::Ref< C, T, F >::isAvailable(), edm::Ref< C, T, F >::isNonnull(), cmsBatch::log, reco::LeafCandidate::p(), reco::LeafCandidate::pt(), reco::GsfElectron::r9(), mps_fire::result, reco::GsfElectron::sigmaEtaEta(), reco::GsfElectron::sigmaIetaIeta(), reco::GsfElectron::sigmaIphiIphi(), reco::GsfElectron::superCluster(), reco::GsfElectron::trackMomentumAtEleClus(), and reco::GsfElectron::trackMomentumAtVtx().

30  {
31  float vars[25];
32 
33  vars[0] = myElectron.fbrem();// fbrem
34  vars[1] = myElectron.eSuperClusterOverP(); //EtotOvePin
35  vars[2] = myElectron.eEleClusterOverPout(); //eleEoPout
36 
37  float etot = myElectron.eSuperClusterOverP()*myElectron.trackMomentumAtVtx().R();
38  float eEcal = myElectron.eEleClusterOverPout()*myElectron.trackMomentumAtEleClus().R();
39  float dP = myElectron.trackMomentumAtVtx().R()-myElectron.trackMomentumAtEleClus().R();
40  vars[3] = (etot-eEcal)/dP; //EBremOverDeltaP
41  vars[4] = std::log(myElectron.sigmaEtaEta()); //logSigmaEtaEta
42  vars[5] = myElectron.deltaEtaEleClusterTrackAtCalo(); //DeltaEtaTrackEcalSeed
43  vars[6] = myElectron.hcalOverEcalBc(); //HoE
44 
45  bool validKF= false;
46  reco::TrackRef myTrackRef = myElectron.closestCtfTrackRef();
47  validKF = (myTrackRef.isAvailable() && myTrackRef.isNonnull());
48  vars[7] = myElectron.gsfTrack()->normalizedChi2(); //gsfchi2
49  vars[8] = (validKF) ? myTrackRef->normalizedChi2() : 0 ; //kfchi2
50  vars[9] = (validKF) ? myTrackRef->hitPattern().trackerLayersWithMeasurement() : -1. ; //kfhits
51 
52  vars[10] = myElectron.gsfTrack().get()->ptModeError()/myElectron.gsfTrack().get()->ptMode() ; //SigmaPtOverPt
53  vars[11] = myElectron.deltaEtaSuperClusterTrackAtVtx(); //deta
54  vars[12] = myElectron.deltaPhiSuperClusterTrackAtVtx(); //dphi
55  vars[13] = myElectron.deltaEtaSeedClusterTrackAtCalo(); //detacalo
56  vars[14] = myElectron.sigmaIetaIeta(); //see
57  vars[15] = myElectron.sigmaIphiIphi(); //spp
58  vars[16] = myElectron.r9(); //R9
59  vars[17] = myElectron.superCluster()->etaWidth(); //etawidth
60  vars[18] = myElectron.superCluster()->phiWidth(); //phiwidth
61  vars[19] = (myElectron.e5x5()) !=0. ? 1.-(myElectron.e1x5()/myElectron.e5x5()) : -1. ; //OneMinusE1x5E5x5
62  vars[20] = (1.0/myElectron.ecalEnergy()) - (1.0 / myElectron.p()); // IoEmIoP
63  vars[21] = myElectron.superCluster()->preshowerEnergy() / myElectron.superCluster()->rawEnergy(); //PreShowerOverRaw
64  vars[22] = pvc.size(); // nPV
65  vars[23] = myElectron.pt(); //pt
66  vars[24] = myElectron.eta(); //eta
67 
68 /*
69  std::cout<<"fbrem "<<fbrem<<std::endl;
70  std::cout<<"EtotOvePin "<<EtotOvePin<<std::endl;
71  std::cout<<"eleEoPout "<<eleEoPout<<std::endl;
72  std::cout<<"EBremOverDeltaP "<<EBremOverDeltaP<<std::endl;
73  std::cout<<"logSigmaEtaEta "<<logSigmaEtaEta<<std::endl;
74  std::cout<<"DeltaEtaTrackEcalSeed "<<DeltaEtaTrackEcalSeed<<std::endl;
75  std::cout<<"HoE "<<HoE<<std::endl;
76  std::cout<<"kfchi2 "<<kfchi2<<std::endl;
77  std::cout<<"kfhits "<<kfhits<<std::endl;
78  std::cout<<"gsfchi2 "<<gsfchi2<<std::endl;
79  std::cout<<"SigmaPtOverPt "<<SigmaPtOverPt<<std::endl;
80  std::cout<<"deta "<<deta<<std::endl;
81  std::cout<<"dphi "<<dphi<<std::endl;
82  std::cout<<"detacalo "<<detacalo<<std::endl;
83  std::cout<<"see "<<see<<std::endl;
84  std::cout<< "spp " << spp<< std::endl;
85  std::cout<< "R9 " << R9<< std::endl;
86  std::cout<< "IoEmIoP " << IoEmIoP<< std::endl;
87  std::cout<<"etawidth "<<etawidth<<std::endl;
88  std::cout<<"phiwidth "<<phiwidth<<std::endl;
89  std::cout<<"OneMinusE1x5E5x5 "<<OneMinusE1x5E5x5<<std::endl;
90  std::cout<<"PreShowerOverRaw "<<PreShowerOverRaw<<std::endl;
91 */
92  bindVariables(vars);
93 
94  double result= gbr_[0]->GetClassifier(vars);
95 
96  return result;
97 }
float sigmaIphiIphi() const
Definition: GsfElectron.h:441
bool isAvailable() const
Definition: Ref.h:575
GsfTrackRef gsfTrack() const override
reference to a GsfTrack
Definition: GsfElectron.h:186
virtual TrackRef closestCtfTrackRef() const
Definition: GsfElectron.h:205
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:251
double eta() const final
momentum pseudorapidity
float eSuperClusterOverP() const
Definition: GsfElectron.h:249
math::XYZVectorF trackMomentumAtVtx() const
Definition: GsfElectron.h:295
double pt() const final
transverse momentum
float fbrem() const
Definition: GsfElectron.h:772
float deltaEtaSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:253
float sigmaIetaIeta() const
Definition: GsfElectron.h:440
float deltaPhiSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:256
float eEleClusterOverPout() const
Definition: GsfElectron.h:252
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:243
float hcalOverEcalBc() const
Definition: GsfElectron.h:452
std::vector< std::unique_ptr< const GBRForest > > gbr_
double p() const final
magnitude of momentum vector
float deltaEtaEleClusterTrackAtCalo() const
Definition: GsfElectron.h:255
float e1x5() const
Definition: GsfElectron.h:442
void bindVariables(float vars[25]) const
float ecalEnergy() const
Definition: GsfElectron.h:859
math::XYZVectorF trackMomentumAtEleClus() const
Definition: GsfElectron.h:298
float e5x5() const
Definition: GsfElectron.h:444
float r9() const
Definition: GsfElectron.h:445
float deltaEtaSeedClusterTrackAtCalo() const
Definition: GsfElectron.h:254
SuperClusterRef superCluster() const override
reference to a SuperCluster
Definition: GsfElectron.h:185
vars
Definition: DeepTauId.cc:77
float sigmaEtaEta() const
Definition: GsfElectron.h:439

Member Data Documentation

const Configuration SoftElectronMVAEstimator::cfg_
private

Definition at line 31 of file SoftElectronMVAEstimator.h.

Referenced by SoftElectronMVAEstimator().

unsigned int SoftElectronMVAEstimator::ExpectedNBins = 1
static

Definition at line 16 of file SoftElectronMVAEstimator.h.

Referenced by SoftElectronMVAEstimator().

std::vector<std::unique_ptr<const GBRForest> > SoftElectronMVAEstimator::gbr_
private

Definition at line 32 of file SoftElectronMVAEstimator.h.

Referenced by mva(), and SoftElectronMVAEstimator().