CMS 3D CMS Logo

Public Member Functions | Private Attributes

PTDRElectronID Class Reference

#include <PTDRElectronID.h>

Inheritance diagram for PTDRElectronID:
ElectronIDAlgo

List of all members.

Public Member Functions

 PTDRElectronID ()
double result (const reco::GsfElectron *, const edm::Event &, const edm::EventSetup &)
void setup (const edm::ParameterSet &conf)
virtual ~PTDRElectronID ()

Private Attributes

std::vector< int > acceptCracks_
edm::ParameterSet cuts_
std::string quality_
std::vector< int > useBremFraction_
std::vector< int > useDeltaEtaIn_
std::vector< int > useDeltaPhiIn_
std::vector< int > useDeltaPhiOut_
std::vector< int > useE9overE25_
std::vector< int > useEoverPIn_
std::vector< int > useEoverPOut_
std::vector< int > useHoverE_
std::vector< int > useInvEMinusInvP_
std::vector< int > useSigmaEtaEta_
std::vector< int > useSigmaPhiPhi_
int variables_

Detailed Description

Definition at line 6 of file PTDRElectronID.h.


Constructor & Destructor Documentation

PTDRElectronID::PTDRElectronID ( ) [inline]

Definition at line 10 of file PTDRElectronID.h.

{};
virtual PTDRElectronID::~PTDRElectronID ( ) [inline, virtual]

Definition at line 12 of file PTDRElectronID.h.

{};

Member Function Documentation

double PTDRElectronID::result ( const reco::GsfElectron electron,
const edm::Event e,
const edm::EventSetup es 
) [virtual]

Reimplemented from ElectronIDAlgo.

Definition at line 37 of file PTDRElectronID.cc.

References acceptCracks_, reco::GsfElectron::BIGBREM, reco::GsfElectron::caloEnergy(), reco::GsfElectron::classification(), cuts_, reco::GsfElectron::deltaEtaSuperClusterTrackAtVtx(), reco::GsfElectron::deltaPhiSeedClusterTrackAtCalo(), reco::GsfElectron::deltaPhiSuperClusterTrackAtVtx(), EcalClusterLazyTools::e3x3(), EcalClusterLazyTools::e5x5(), reco::GsfElectron::eSeedClusterOverPout(), reco::GsfElectron::eSuperClusterOverP(), reco::GsfElectron::GAP, ElectronIDAlgo::getClusterShape(), edm::ParameterSet::getParameter(), reco::GsfElectron::GOLDEN, reco::GsfElectron::hadronicOverEm(), reco::GsfElectron::isEB(), reco::GsfElectron::isEE(), EcalClusterLazyTools::localCovariances(), reco::GsfElectron::SHOWERING, mathSSE::sqrt(), reco::GsfElectron::superCluster(), reco::GsfElectron::trackMomentumAtVtx(), reco::GsfElectron::trackMomentumOut(), reco::GsfElectron::UNKNOWN, useBremFraction_, useDeltaEtaIn_, useDeltaPhiIn_, useDeltaPhiOut_, useE9overE25_, useEoverPIn_, useEoverPOut_, useHoverE_, useInvEMinusInvP_, useSigmaEtaEta_, useSigmaPhiPhi_, relativeConstraints::value, and variables_.

                                                       {

  //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::NARROW)    icut=2;
       if (elClass == reco::GsfElectron::SHOWERING) icut=3;
       if (elClass == reco::GsfElectron::GAP)       icut=8;
     }
  if (electron->isEE()) //endcap
     {
       if (elClass == reco::GsfElectron::GOLDEN)    icut=4;
       if (elClass == reco::GsfElectron::BIGBREM)   icut=5;
       //if (elClass == reco::GsfElectron::NARROW)    icut=6;
       if (elClass == reco::GsfElectron::SHOWERING) icut=7;
       if (elClass == reco::GsfElectron::GAP)       icut=8;
     }
  if (elClass == reco::GsfElectron::UNKNOWN) 
     {
       edm::LogError("PTDRElectronID") << "Error: unrecognized electron classification ";
       return 1.;
     }

  if (acceptCracks_[variables_])
    if (elClass == reco::GsfElectron::GAP) return 1.;
  
  if (useEoverPIn_[variables_]) {
    double value = electron->eSuperClusterOverP();
    std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("EoverPInMax");
    std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("EoverPInMin");
    if (value<mincut[icut] || value>maxcut[icut]) return 0.;
  }

  if (useDeltaEtaIn_[variables_]) {
    double value = electron->deltaEtaSuperClusterTrackAtVtx();
    std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("deltaEtaIn");
    if (fabs(value)>maxcut[icut]) return 0.;
  }

  if (useDeltaPhiIn_[variables_]) {
    double value = electron->deltaPhiSuperClusterTrackAtVtx();
    std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("deltaPhiIn");
    if (fabs(value)>maxcut[icut]) return 0.;
  }

  if (useHoverE_[variables_]) {
    double value = electron->hadronicOverEm();
    std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("HoverE");
    if (value>maxcut[icut]) return 0.;
  }

  if (useEoverPOut_[variables_]) {
    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.;
  }

  if (useDeltaPhiOut_[variables_]) {
    double value = electron->deltaPhiSeedClusterTrackAtCalo();
    std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("deltaPhiOut");
    if (fabs(value)>maxcut[icut]) return 0.;
  }

  if (useInvEMinusInvP_[variables_]) {
    double value = (1./electron->caloEnergy())-(1./electron->trackMomentumAtVtx().R());
    std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("invEMinusInvP");
    if (value>maxcut[icut]) return 0.;
  }

  if (useBremFraction_[variables_]) {
    double value = electron->trackMomentumAtVtx().R()-electron->trackMomentumOut().R();
    std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("bremFraction");
    if (value<mincut[icut]) return 0.;
  }

  EcalClusterLazyTools lazyTools = getClusterShape(e,es);
  std::vector<float> vCov = lazyTools.localCovariances(*(electron->superCluster()->seed())) ;
  //std::vector<float> vCov = lazyTools.covariances(*(electron->superCluster()->seed())) ;
    
  if (useE9overE25_[variables_]) {
    double value = lazyTools.e3x3(*(electron->superCluster()->seed()))/lazyTools.e5x5(*(electron->superCluster()->seed()));
    std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("E9overE25");
    if (fabs(value)<mincut[icut]) return 0.;
  }

  if (useSigmaEtaEta_[variables_]) {
    std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("sigmaEtaEtaMax");
    std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("sigmaEtaEtaMin");
    if (sqrt(vCov[0])<mincut[icut] || sqrt(vCov[0])>maxcut[icut]) return 0.;
  }

  if (useSigmaPhiPhi_[variables_]) {
    std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("sigmaPhiPhiMin");
    std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("sigmaPhiPhiMax");
    if (sqrt(vCov[1])<mincut[icut] || sqrt(vCov[1])>maxcut[icut]) return 0.;
  }

  return 1.;

}
void PTDRElectronID::setup ( const edm::ParameterSet conf) [virtual]

Reimplemented from ElectronIDAlgo.

Definition at line 3 of file PTDRElectronID.cc.

References acceptCracks_, ElectronIDAlgo::baseSetup(), cuts_, Exception, edm::ParameterSet::getParameter(), quality_, useBremFraction_, useDeltaEtaIn_, useDeltaPhiIn_, useDeltaPhiOut_, useE9overE25_, useEoverPIn_, useEoverPOut_, useHoverE_, useInvEMinusInvP_, useSigmaEtaEta_, useSigmaPhiPhi_, and variables_.

                                                      {

  // Get all the parameters
  baseSetup(conf);
  
  quality_ =  conf.getParameter<std::string>("electronQuality");
  
  useEoverPIn_ = conf.getParameter<std::vector<int> >("useEoverPIn");
  useDeltaEtaIn_ = conf.getParameter<std::vector<int> >("useDeltaEtaIn");
  useDeltaPhiIn_ = conf.getParameter<std::vector<int> >("useDeltaPhiIn");
  useHoverE_ = conf.getParameter<std::vector<int> >("useHoverE");
  useE9overE25_ = conf.getParameter<std::vector<int> >("useE9overE25");
  useEoverPOut_ = conf.getParameter<std::vector<int> >("useEoverPOut");
  useDeltaPhiOut_ = conf.getParameter<std::vector<int> >("useDeltaPhiOut");
  useInvEMinusInvP_ = conf.getParameter<std::vector<int> >("useInvEMinusInvP");
  useBremFraction_ = conf.getParameter<std::vector<int> >("useBremFraction");
  useSigmaEtaEta_ = conf.getParameter<std::vector<int> >("useSigmaEtaEta");
  useSigmaPhiPhi_ = conf.getParameter<std::vector<int> >("useSigmaPhiPhi");
  acceptCracks_ = conf.getParameter<std::vector<int> >("acceptCracks");
  
  if (quality_=="tight") {
    cuts_ = conf.getParameter<edm::ParameterSet>("tightEleIDCuts");
    variables_ = 2 ;
  } else if (quality_=="medium") {
    cuts_ = conf.getParameter<edm::ParameterSet>("mediumEleIDCuts");
    variables_ = 1 ;
  } else if (quality_=="loose") {
    cuts_ = conf.getParameter<edm::ParameterSet>("looseEleIDCuts");
    variables_ = 0 ;
  } else {
     throw cms::Exception("Configuration") << "Invalid electronQuality parameter in PTDElectronID: must be tight, medium or loose." ;
  }  
}

Member Data Documentation

std::vector<int> PTDRElectronID::acceptCracks_ [private]

Definition at line 32 of file PTDRElectronID.h.

Referenced by result(), and setup().

Definition at line 34 of file PTDRElectronID.h.

Referenced by result(), and setup().

std::string PTDRElectronID::quality_ [private]

Definition at line 19 of file PTDRElectronID.h.

Referenced by setup().

std::vector<int> PTDRElectronID::useBremFraction_ [private]

Definition at line 29 of file PTDRElectronID.h.

Referenced by result(), and setup().

std::vector<int> PTDRElectronID::useDeltaEtaIn_ [private]

Definition at line 22 of file PTDRElectronID.h.

Referenced by result(), and setup().

std::vector<int> PTDRElectronID::useDeltaPhiIn_ [private]

Definition at line 23 of file PTDRElectronID.h.

Referenced by result(), and setup().

std::vector<int> PTDRElectronID::useDeltaPhiOut_ [private]

Definition at line 27 of file PTDRElectronID.h.

Referenced by result(), and setup().

std::vector<int> PTDRElectronID::useE9overE25_ [private]

Definition at line 25 of file PTDRElectronID.h.

Referenced by result(), and setup().

std::vector<int> PTDRElectronID::useEoverPIn_ [private]

Definition at line 21 of file PTDRElectronID.h.

Referenced by result(), and setup().

std::vector<int> PTDRElectronID::useEoverPOut_ [private]

Definition at line 26 of file PTDRElectronID.h.

Referenced by result(), and setup().

std::vector<int> PTDRElectronID::useHoverE_ [private]

Definition at line 24 of file PTDRElectronID.h.

Referenced by result(), and setup().

std::vector<int> PTDRElectronID::useInvEMinusInvP_ [private]

Definition at line 28 of file PTDRElectronID.h.

Referenced by result(), and setup().

std::vector<int> PTDRElectronID::useSigmaEtaEta_ [private]

Definition at line 30 of file PTDRElectronID.h.

Referenced by result(), and setup().

std::vector<int> PTDRElectronID::useSigmaPhiPhi_ [private]

Definition at line 31 of file PTDRElectronID.h.

Referenced by result(), and setup().

Definition at line 36 of file PTDRElectronID.h.

Referenced by result(), and setup().