#include <ElectronMVAEstimator.h>
Public Member Functions | |
ElectronMVAEstimator () | |
ElectronMVAEstimator (std::string fileName) | |
double | mva (const reco::GsfElectron &myElectron, int nvertices=0) |
~ElectronMVAEstimator () | |
Private Member Functions | |
void | bindVariables () |
void | init (std::string fileName) |
Private Attributes | |
Float_t | absdcot |
Float_t | absdist |
Float_t | dcot |
Float_t | detaeleout |
Float_t | detain |
Float_t | dist |
Float_t | dphiin |
Float_t | e1x5e5x5 |
Int_t | ecalseed |
Float_t | eleopout |
Float_t | eop |
Float_t | eta |
Float_t | fbrem |
Float_t | hoe |
Float_t | kfchi2 |
Int_t | kfhits |
Int_t | mishits |
Float_t | mykfhits |
Float_t | mymishits |
Float_t | myNvtx |
Int_t | Nvtx |
Float_t | pt |
Float_t | sieie |
TMVA::Reader * | tmvaReader_ |
Definition at line 8 of file ElectronMVAEstimator.h.
ElectronMVAEstimator::ElectronMVAEstimator | ( | ) |
Definition at line 8 of file ElectronMVAEstimator.cc.
{}
ElectronMVAEstimator::ElectronMVAEstimator | ( | std::string | fileName | ) |
ElectronMVAEstimator::~ElectronMVAEstimator | ( | ) | [inline] |
Definition at line 12 of file ElectronMVAEstimator.h.
{;}
void ElectronMVAEstimator::bindVariables | ( | ) | [private] |
Definition at line 101 of file ElectronMVAEstimator.cc.
References absdcot, absdist, dcot, detaeleout, detain, dist, dphiin, e1x5e5x5, eleopout, eop, fbrem, kfchi2, kfhits, mishits, mykfhits, mymishits, myNvtx, and Nvtx.
Referenced by mva().
{ if(fbrem < -1.) fbrem = -1.; detain = fabs(detain); if(detain > 0.06) detain = 0.06; dphiin = fabs(dphiin); if(dphiin > 0.6) dphiin = 0.6; if(eop > 20.) eop = 20.; if(eleopout > 20.) eleopout = 20; detaeleout = fabs(detaeleout); if(detaeleout > 0.2) detaeleout = 0.2; mykfhits = float(kfhits); mymishits = float(mishits); if(kfchi2 < 0.) kfchi2 = 0.; if(kfchi2 > 15.) kfchi2 = 15.; if(e1x5e5x5 < -1.) e1x5e5x5 = -1; if(e1x5e5x5 > 2.) e1x5e5x5 = 2.; if(dist > 15.) dist = 15.; if(dist < -15.) dist = -15.; if(dcot > 3.) dcot = 3.; if(dcot < -3.) dcot = -3.; absdist = fabs(dist); absdcot = fabs(dcot); myNvtx = float(Nvtx); }
void ElectronMVAEstimator::init | ( | std::string | fileName | ) | [private] |
Definition at line 15 of file ElectronMVAEstimator.cc.
References absdcot, absdist, detaeleout, detain, dphiin, e1x5e5x5, ecalseed, eleopout, eop, eta, fbrem, hoe, kfchi2, mykfhits, mymishits, myNvtx, pt, sieie, and tmvaReader_.
Referenced by ElectronMVAEstimator().
{ tmvaReader_ = new TMVA::Reader("!Color:Silent"); tmvaReader_->AddVariable("fbrem",&fbrem); tmvaReader_->AddVariable("detain", &detain); tmvaReader_->AddVariable("dphiin", &dphiin); tmvaReader_->AddVariable("sieie", &sieie); tmvaReader_->AddVariable("hoe", &hoe); tmvaReader_->AddVariable("eop", &eop); tmvaReader_->AddVariable("e1x5e5x5", &e1x5e5x5); tmvaReader_->AddVariable("eleopout", &eleopout); tmvaReader_->AddVariable("detaeleout", &detaeleout); tmvaReader_->AddVariable("kfchi2", &kfchi2); tmvaReader_->AddVariable("kfhits", &mykfhits); tmvaReader_->AddVariable("mishits",&mymishits); tmvaReader_->AddVariable("dist", &absdist); tmvaReader_->AddVariable("dcot", &absdcot); tmvaReader_->AddVariable("nvtx", &myNvtx); tmvaReader_->AddSpectator("eta",&eta); tmvaReader_->AddSpectator("pt",&pt); tmvaReader_->AddSpectator("ecalseed",&ecalseed); // Taken from Daniele (his mail from the 30/11) // tmvaReader_->BookMVA("BDTSimpleCat","../Training/weights_Root527b_3Depth_DanVarConvRej_2PtBins_10Pt_800TPrune5_Min100Events_NoBjets_half/TMVA_BDTSimpleCat.weights.xm"); // training of the 7/12 with Nvtx added tmvaReader_->BookMVA("BDTSimpleCat",fileName.c_str()); }
double ElectronMVAEstimator::mva | ( | const reco::GsfElectron & | myElectron, |
int | nvertices = 0 |
||
) |
Definition at line 45 of file ElectronMVAEstimator.cc.
References bindVariables(), reco::GsfElectron::closestCtfTrackRef(), reco::GsfElectron::convDcot(), reco::GsfElectron::convDist(), dcot, reco::GsfElectron::deltaEtaEleClusterTrackAtCalo(), reco::GsfElectron::deltaEtaSuperClusterTrackAtVtx(), reco::GsfElectron::deltaPhiSuperClusterTrackAtVtx(), detaeleout, detain, dist, dphiin, reco::GsfElectron::e1x5(), e1x5e5x5, reco::GsfElectron::e5x5(), reco::GsfElectron::ecalDrivenSeed(), ecalseed, reco::GsfElectron::eEleClusterOverPout(), eleopout, eop, reco::GsfElectron::eSuperClusterOverP(), reco::LeafCandidate::eta(), eta, reco::GsfElectron::fbrem(), fbrem, reco::GsfElectron::gsfTrack(), reco::GsfElectron::hcalOverEcal(), hoe, edm::Ref< C, T, F >::isAvailable(), edm::Ref< C, T, F >::isNonnull(), kfchi2, kfhits, mishits, Nvtx, pt, reco::LeafCandidate::pt(), query::result, sieie, reco::GsfElectron::sigmaIetaIeta(), and tmvaReader_.
{ fbrem = myElectron.fbrem(); detain = myElectron.deltaEtaSuperClusterTrackAtVtx(); dphiin = myElectron.deltaPhiSuperClusterTrackAtVtx(); sieie = myElectron.sigmaIetaIeta(); hoe = myElectron.hcalOverEcal(); eop = myElectron.eSuperClusterOverP(); e1x5e5x5 = (myElectron.e5x5()) !=0. ? 1.-(myElectron.e1x5()/myElectron.e5x5()) : -1. ; eleopout = myElectron.eEleClusterOverPout(); detaeleout = myElectron.deltaEtaEleClusterTrackAtCalo(); bool validKF= false; reco::TrackRef myTrackRef = myElectron.closestCtfTrackRef(); validKF = (myTrackRef.isAvailable()); validKF = (myTrackRef.isNonnull()); kfchi2 = (validKF) ? myTrackRef->normalizedChi2() : 0 ; kfhits = (validKF) ? myTrackRef->hitPattern().trackerLayersWithMeasurement() : -1. ; dist = myElectron.convDist(); dcot = myElectron.convDcot(); eta = myElectron.eta(); pt = myElectron.pt(); mishits = myElectron.gsfTrack()->trackerExpectedHitsInner().numberOfLostHits(); ecalseed = myElectron.ecalDrivenSeed(); Nvtx = nvertices; bindVariables(); double result = tmvaReader_->EvaluateMVA("BDTSimpleCat"); // // std::cout << "fbrem" <<fbrem << std::endl; // std::cout << "detain"<< detain << std::endl; // std::cout << "dphiin"<< dphiin << std::endl; // std::cout << "sieie"<< sieie << std::endl; // std::cout << "hoe"<< hoe << std::endl; // std::cout << "eop"<< eop << std::endl; // std::cout << "e1x5e5x5"<< e1x5e5x5 << std::endl; // std::cout << "eleopout"<< eleopout << std::endl; // std::cout << "detaeleout"<< detaeleout << std::endl; // std::cout << "kfchi2"<< kfchi2 << std::endl; // std::cout << "kfhits"<< mykfhits << std::endl; // std::cout << "mishits"<<mymishits << std::endl; // std::cout << "dist"<< absdist << std::endl; // std::cout << "dcot"<< absdcot << std::endl; // // std::cout << "eta"<<eta << std::endl; // std::cout << "pt"<<pt << std::endl; // std::cout << "ecalseed"<<ecalseed << std::endl; // // std::cout << " MVA " << result << std::endl; return result; }
Float_t ElectronMVAEstimator::absdcot [private] |
Definition at line 42 of file ElectronMVAEstimator.h.
Referenced by bindVariables(), and init().
Float_t ElectronMVAEstimator::absdist [private] |
Definition at line 41 of file ElectronMVAEstimator.h.
Referenced by bindVariables(), and init().
Float_t ElectronMVAEstimator::dcot [private] |
Definition at line 33 of file ElectronMVAEstimator.h.
Referenced by bindVariables(), and mva().
Float_t ElectronMVAEstimator::detaeleout [private] |
Definition at line 30 of file ElectronMVAEstimator.h.
Referenced by bindVariables(), init(), and mva().
Float_t ElectronMVAEstimator::detain [private] |
Definition at line 23 of file ElectronMVAEstimator.h.
Referenced by bindVariables(), init(), and mva().
Float_t ElectronMVAEstimator::dist [private] |
Definition at line 32 of file ElectronMVAEstimator.h.
Referenced by bindVariables(), and mva().
Float_t ElectronMVAEstimator::dphiin [private] |
Definition at line 24 of file ElectronMVAEstimator.h.
Referenced by bindVariables(), init(), and mva().
Float_t ElectronMVAEstimator::e1x5e5x5 [private] |
Definition at line 28 of file ElectronMVAEstimator.h.
Referenced by bindVariables(), init(), and mva().
Int_t ElectronMVAEstimator::ecalseed [private] |
Definition at line 38 of file ElectronMVAEstimator.h.
Float_t ElectronMVAEstimator::eleopout [private] |
Definition at line 29 of file ElectronMVAEstimator.h.
Referenced by bindVariables(), init(), and mva().
Float_t ElectronMVAEstimator::eop [private] |
Definition at line 27 of file ElectronMVAEstimator.h.
Referenced by bindVariables(), init(), and mva().
Float_t ElectronMVAEstimator::eta [private] |
Definition at line 34 of file ElectronMVAEstimator.h.
Float_t ElectronMVAEstimator::fbrem [private] |
Definition at line 22 of file ElectronMVAEstimator.h.
Referenced by bindVariables(), init(), and mva().
Float_t ElectronMVAEstimator::hoe [private] |
Definition at line 26 of file ElectronMVAEstimator.h.
Float_t ElectronMVAEstimator::kfchi2 [private] |
Definition at line 31 of file ElectronMVAEstimator.h.
Referenced by bindVariables(), init(), and mva().
Int_t ElectronMVAEstimator::kfhits [private] |
Definition at line 36 of file ElectronMVAEstimator.h.
Referenced by bindVariables(), and mva().
Int_t ElectronMVAEstimator::mishits [private] |
Definition at line 37 of file ElectronMVAEstimator.h.
Referenced by bindVariables(), and mva().
Float_t ElectronMVAEstimator::mykfhits [private] |
Definition at line 43 of file ElectronMVAEstimator.h.
Referenced by bindVariables(), and init().
Float_t ElectronMVAEstimator::mymishits [private] |
Definition at line 44 of file ElectronMVAEstimator.h.
Referenced by bindVariables(), and init().
Float_t ElectronMVAEstimator::myNvtx [private] |
Definition at line 45 of file ElectronMVAEstimator.h.
Referenced by bindVariables(), and init().
Int_t ElectronMVAEstimator::Nvtx [private] |
Definition at line 39 of file ElectronMVAEstimator.h.
Referenced by bindVariables(), and mva().
Float_t ElectronMVAEstimator::pt [private] |
Definition at line 35 of file ElectronMVAEstimator.h.
Float_t ElectronMVAEstimator::sieie [private] |
Definition at line 25 of file ElectronMVAEstimator.h.
TMVA::Reader* ElectronMVAEstimator::tmvaReader_ [private] |
Definition at line 20 of file ElectronMVAEstimator.h.