#include <ClassBasedElectronID.h>
Public Member Functions | |
ClassBasedElectronID () | |
double | result (const reco::GsfElectron *, const edm::Event &, const edm::EventSetup &) |
void | setup (const edm::ParameterSet &conf) |
virtual | ~ClassBasedElectronID () |
Private Attributes | |
edm::ParameterSet | cuts_ |
std::string | quality_ |
Definition at line 6 of file ClassBasedElectronID.h.
ClassBasedElectronID::ClassBasedElectronID | ( | ) | [inline] |
Definition at line 10 of file ClassBasedElectronID.h.
{};
virtual ClassBasedElectronID::~ClassBasedElectronID | ( | ) | [inline, virtual] |
Definition at line 12 of file ClassBasedElectronID.h.
{};
double ClassBasedElectronID::result | ( | const reco::GsfElectron * | electron, |
const edm::Event & | e, | ||
const edm::EventSetup & | es | ||
) | [virtual] |
Reimplemented from ElectronIDAlgo.
Definition at line 29 of file ClassBasedElectronID.cc.
References reco::GsfElectron::BIGBREM, reco::LeafCandidate::charge(), reco::GsfElectron::classification(), cuts_, reco::GsfElectron::deltaEtaSuperClusterTrackAtVtx(), reco::GsfElectron::deltaPhiSuperClusterTrackAtVtx(), reco::GsfElectron::eSeedClusterOverPout(), reco::GsfElectron::GAP, edm::ParameterSet::getParameter(), reco::GsfElectron::GOLDEN, reco::GsfElectron::hadronicOverEm(), reco::GsfElectron::isEB(), reco::GsfElectron::isEE(), reco::GsfElectron::SHOWERING, reco::GsfElectron::sigmaIetaIeta(), reco::GsfElectron::UNKNOWN, and relativeConstraints::value.
{ //determine which element of the cut arrays in cfi file to read //depending on the electron classification int icut=0; int elClass = electron->classification() ; if (electron->isEB()) //barrel { if (elClass == reco::GsfElectron::GOLDEN) icut=0; if (elClass == reco::GsfElectron::BIGBREM) icut=1; if (elClass == reco::GsfElectron::SHOWERING) icut=2; if (elClass == reco::GsfElectron::GAP) icut=6; } if (electron->isEE()) //endcap { if (elClass == reco::GsfElectron::GOLDEN) icut=3; if (elClass == reco::GsfElectron::BIGBREM) icut=4; if (elClass == reco::GsfElectron::SHOWERING) icut=5; if (elClass == reco::GsfElectron::GAP) icut=7; } if (elClass == reco::GsfElectron::UNKNOWN) { edm::LogError("ClassBasedElectronID") << "Error: unrecognized electron classification "; return 1.; } bool useDeltaEtaIn = true; bool useSigmaIetaIeta = true; bool useHoverE = true; bool useEoverPOut = true; bool useDeltaPhiInCharge = true; // DeltaEtaIn if (useDeltaEtaIn) { double value = electron->deltaEtaSuperClusterTrackAtVtx(); std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("deltaEtaIn"); if (fabs(value)>maxcut[icut]) return 0.; } // SigmaIetaIeta if(useSigmaIetaIeta) { double value = electron->sigmaIetaIeta() ; std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("sigmaIetaIetaMax"); std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("sigmaIetaIetaMin"); if(value<mincut[icut] || value>maxcut[icut]) return 0.; } // H/E if (useHoverE) { //_[variables_]) { double value = electron->hadronicOverEm(); std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("HoverE"); if (value>maxcut[icut]) return 0.; } // if use // Eseed/Pout if (useEoverPOut) { double value = electron->eSeedClusterOverPout(); std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("EoverPOutMax"); std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("EoverPOutMin"); if (value<mincut[icut] || value>maxcut[icut]) return 0.; } // DeltaPhiIn*Charge if (useDeltaPhiInCharge) { double value1 = electron->deltaPhiSuperClusterTrackAtVtx(); double value2 = electron->charge(); double value = value1*value2; std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("deltaPhiInChargeMax"); std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("deltaPhiInChargeMin"); if (value<mincut[icut] || value>maxcut[icut]) return 0.; } return 1.; }
void ClassBasedElectronID::setup | ( | const edm::ParameterSet & | conf | ) | [virtual] |
Reimplemented from ElectronIDAlgo.
Definition at line 4 of file ClassBasedElectronID.cc.
References ElectronIDAlgo::baseSetup(), cuts_, cmsRelvalreport::exit, edm::ParameterSet::getParameter(), and quality_.
{ // Get all the parameters baseSetup(conf); quality_ = conf.getParameter<std::string>("electronQuality"); if(quality_=="Eff95Cuts") { cuts_ = conf.getParameter<edm::ParameterSet>("Eff95Cuts"); } else if(quality_=="Eff90Cuts") { cuts_ = conf.getParameter<edm::ParameterSet>("Eff90Cuts"); } else { edm::LogError("ClassBasedElectronID") << "Invalid electronQuality parameter: must be tight, medium or loose." ; exit (1); } } // end of setup
edm::ParameterSet ClassBasedElectronID::cuts_ [private] |
Definition at line 34 of file ClassBasedElectronID.h.
std::string ClassBasedElectronID::quality_ [private] |
Definition at line 19 of file ClassBasedElectronID.h.
Referenced by setup().