CMS 3D CMS Logo

List of all members | Public Member Functions | Public Attributes
ZElectronsSelector Class Reference

Public Member Functions

const float getEffectiveArea (float eta) const
 
void newEvent (const edm::Event &, const edm::EventSetup &)
 
bool operator() (const reco::GsfElectron &) const
 
void printEffectiveAreas () const
 
 ZElectronsSelector (const edm::ParameterSet &, edm::ConsumesCollector &iC)
 

Public Attributes

edm::Handle< double > _rhoHandle
 
std::vector< double > absEtaMax_
 
std::vector< double > absEtaMin_
 
vector< double > dEtaInSeedCut
 
vector< double > dPhiInCut
 
std::vector< double > effectiveAreaValues_
 
vector< double > EInvMinusPInv
 
edm::ParameterSet eleIDWP
 
vector< double > hOverECut
 
vector< int > missHits
 
vector< double > relCombIso
 
vector< double > sigmaIEtaIEtaCut
 
edm::EDGetTokenT< reco::GsfElectronCollectiontheGsfEToken
 
edm::EDGetTokenT< double > theRhoToken
 

Detailed Description

Definition at line 45 of file ZElectronsSelector.cc.

Constructor & Destructor Documentation

ZElectronsSelector::ZElectronsSelector ( const edm::ParameterSet cfg,
edm::ConsumesCollector iC 
)

Definition at line 92 of file ZElectronsSelector.cc.

References absEtaMax_, absEtaMin_, dEtaInSeedCut, dPhiInCut, effectiveAreaValues_, EInvMinusPInv, eleIDWP, edm::ParameterSet::getParameter(), hOverECut, missHits, relCombIso, and sigmaIEtaIEtaCut.

93  : theRhoToken(iC.consumes<double>(cfg.getParameter<edm::InputTag>("rho"))) {
94  absEtaMin_ = cfg.getParameter<std::vector<double> >("absEtaMin");
95  absEtaMax_ = cfg.getParameter<std::vector<double> >("absEtaMax");
96  effectiveAreaValues_ = cfg.getParameter<std::vector<double> >("effectiveAreaValues");
97  //printEffectiveAreas();
98  eleIDWP = cfg.getParameter<edm::ParameterSet>("eleID");
99 
100  missHits = eleIDWP.getParameter<std::vector<int> >("missingHitsCut");
101  sigmaIEtaIEtaCut = eleIDWP.getParameter<std::vector<double> >("full5x5_sigmaIEtaIEtaCut");
102  dEtaInSeedCut = eleIDWP.getParameter<std::vector<double> >("dEtaInSeedCut");
103  dPhiInCut = eleIDWP.getParameter<std::vector<double> >("dPhiInCut");
104  hOverECut = eleIDWP.getParameter<std::vector<double> >("hOverECut");
105  relCombIso = eleIDWP.getParameter<std::vector<double> >("relCombIsolationWithEACut");
106  EInvMinusPInv = eleIDWP.getParameter<std::vector<double> >("EInverseMinusPInverseCut");
107 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
vector< double > EInvMinusPInv
std::vector< double > absEtaMax_
vector< double > dEtaInSeedCut
std::vector< double > effectiveAreaValues_
edm::EDGetTokenT< double > theRhoToken
vector< double > relCombIso
edm::ParameterSet eleIDWP
vector< double > dPhiInCut
std::vector< double > absEtaMin_
vector< double > hOverECut
vector< double > sigmaIEtaIEtaCut

Member Function Documentation

const float ZElectronsSelector::getEffectiveArea ( float  eta) const

Definition at line 79 of file ZElectronsSelector.cc.

References funct::abs(), ticl::constants::nEtaBins, and parallelization::uint.

Referenced by operator()().

79  {
80  float effArea = 0;
81  uint nEtaBins = absEtaMin_.size();
82  for (uint iEta = 0; iEta < nEtaBins; iEta++) {
83  if (std::abs(eta) >= absEtaMin_[iEta] && std::abs(eta) < absEtaMax_[iEta]) {
84  effArea = effectiveAreaValues_[iEta];
85  break;
86  }
87  }
88 
89  return effArea;
90 }
std::vector< double > absEtaMax_
std::vector< double > effectiveAreaValues_
constexpr int nEtaBins
Definition: Common.h:11
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< double > absEtaMin_
void ZElectronsSelector::newEvent ( const edm::Event ev,
const edm::EventSetup  
)

Definition at line 109 of file ZElectronsSelector.cc.

References _rhoHandle, edm::Event::getByToken(), and theRhoToken.

109  {
111 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
edm::Handle< double > _rhoHandle
edm::EDGetTokenT< double > theRhoToken
bool ZElectronsSelector::operator() ( const reco::GsfElectron el) const

Definition at line 113 of file ZElectronsSelector.cc.

References _rhoHandle, funct::abs(), reco::GsfElectron::deltaEtaSeedClusterTrackAtVtx(), reco::GsfElectron::deltaPhiSuperClusterTrackAtVtx(), dEtaInSeedCut, dPhiInCut, reco::GsfElectron::ecalEnergy(), EInvMinusPInv, reco::GsfElectron::eSuperClusterOverP(), dqmMemoryStats::float, reco::GsfElectron::full5x5_sigmaIetaIeta(), edm::Ref< C, T, F >::get(), getEffectiveArea(), reco::GsfElectron::gsfTrack(), reco::GsfElectron::hadronicOverEm(), hOverECut, reco::GsfElectron::isEB(), reco::GsfElectron::isEE(), edm::HandleBase::isValid(), SiStripPI::max, missHits, reco::HitPattern::MISSING_INNER_HITS, reco::GsfElectron::pfIsolationVariables(), position, edm::Handle< T >::product(), reco::LeafCandidate::pt(), relCombIso, rho, sigmaIEtaIEtaCut, reco::GsfElectron::PflowIsolationVariables::sumChargedHadronPt, reco::GsfElectron::PflowIsolationVariables::sumNeutralHadronEt, reco::GsfElectron::PflowIsolationVariables::sumPhotonEt, and reco::GsfElectron::superCluster().

113  {
114  float pt_e = el.pt();
115  unsigned int ind = 0;
116  auto etrack = el.gsfTrack();
117  float abseta = fabs((el.superCluster().get())->position().eta());
118 
119  if (el.isEB()) {
120  if (abseta > 1.479)
121  return false; // check if it is really needed
122  }
123  if (el.isEE()) {
124  ind = 1;
125  if (abseta < 1.479)
126  return false; // check if it is really needed
127  if (abseta >= 2.5)
128  return false; // check if it is really needed
129  }
130 
131  if (etrack->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS) > missHits[ind])
132  return false;
133  if (el.full5x5_sigmaIetaIeta() > sigmaIEtaIEtaCut[ind])
134  return false;
135  if (fabs(el.deltaPhiSuperClusterTrackAtVtx()) > dPhiInCut[ind])
136  return false;
137  if (fabs(el.deltaEtaSeedClusterTrackAtVtx()) > dEtaInSeedCut[ind])
138  return false;
139  if (el.hadronicOverEm() > hOverECut[ind])
140  return false;
141  const float eA = getEffectiveArea(abseta);
142  const float rho = _rhoHandle.isValid() ? (float)(*_rhoHandle.product()) : 0;
144  std::max(float(0.0),
146  relCombIso[ind] * pt_e)
147  return false;
148  const float ecal_energy_inverse = 1.0 / el.ecalEnergy();
149  const float eSCoverP = el.eSuperClusterOverP();
150  if (std::abs(1.0 - eSCoverP) * ecal_energy_inverse > EInvMinusPInv[ind])
151  return false;
152 
153  return true;
154 }
const PflowIsolationVariables & pfIsolationVariables() const
Definition: GsfElectron.h:650
GsfTrackRef gsfTrack() const override
reference to a GsfTrack
Definition: GsfElectron.h:156
float eSuperClusterOverP() const
Definition: GsfElectron.h:221
vector< double > EInvMinusPInv
edm::Handle< double > _rhoHandle
vector< double > dEtaInSeedCut
double pt() const final
transverse momentum
bool isEE() const
Definition: GsfElectron.h:329
bool isEB() const
Definition: GsfElectron.h:328
float full5x5_sigmaIetaIeta() const
Definition: GsfElectron.h:435
float hadronicOverEm() const
Definition: GsfElectron.h:468
float sumPhotonEt
sum pt of PF photons // old float photonIso ;
Definition: GsfElectron.h:602
float deltaPhiSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:228
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float sumNeutralHadronEt
sum pt of neutral hadrons // old float neutralHadronIso ;
Definition: GsfElectron.h:601
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:232
bool isValid() const
Definition: HandleBase.h:70
float deltaEtaSeedClusterTrackAtVtx() const
Definition: GsfElectron.h:231
vector< double > relCombIso
T const * product() const
Definition: Handle.h:69
float ecalEnergy() const
Definition: GsfElectron.h:812
static int position[264][3]
Definition: ReadPGInfo.cc:289
SuperClusterRef superCluster() const override
reference to a SuperCluster
Definition: GsfElectron.h:155
vector< double > dPhiInCut
float sumChargedHadronPt
sum-pt of charged Hadron // old float chargedHadronIso ;
Definition: GsfElectron.h:600
const float getEffectiveArea(float eta) const
vector< double > hOverECut
vector< double > sigmaIEtaIEtaCut
void ZElectronsSelector::printEffectiveAreas ( ) const

Definition at line 72 of file ZElectronsSelector.cc.

References ticl::constants::nEtaBins, and parallelization::uint.

72  {
73  printf(" eta_min eta_max effective area\n");
74  uint nEtaBins = absEtaMin_.size();
75  for (uint iEta = 0; iEta < nEtaBins; iEta++) {
76  printf(" %8.4f %8.4f %8.5f\n", absEtaMin_[iEta], absEtaMax_[iEta], effectiveAreaValues_[iEta]);
77  }
78 }
std::vector< double > absEtaMax_
std::vector< double > effectiveAreaValues_
constexpr int nEtaBins
Definition: Common.h:11
std::vector< double > absEtaMin_

Member Data Documentation

edm::Handle<double> ZElectronsSelector::_rhoHandle

Definition at line 55 of file ZElectronsSelector.cc.

Referenced by newEvent(), and operator()().

std::vector<double> ZElectronsSelector::absEtaMax_

Definition at line 58 of file ZElectronsSelector.cc.

Referenced by ZElectronsSelector().

std::vector<double> ZElectronsSelector::absEtaMin_

Definition at line 57 of file ZElectronsSelector.cc.

Referenced by ZElectronsSelector().

vector<double> ZElectronsSelector::dEtaInSeedCut

Definition at line 65 of file ZElectronsSelector.cc.

Referenced by operator()(), and ZElectronsSelector().

vector<double> ZElectronsSelector::dPhiInCut

Definition at line 66 of file ZElectronsSelector.cc.

Referenced by operator()(), and ZElectronsSelector().

std::vector<double> ZElectronsSelector::effectiveAreaValues_

Definition at line 59 of file ZElectronsSelector.cc.

Referenced by ZElectronsSelector().

vector<double> ZElectronsSelector::EInvMinusPInv

Definition at line 69 of file ZElectronsSelector.cc.

Referenced by operator()(), and ZElectronsSelector().

edm::ParameterSet ZElectronsSelector::eleIDWP

Definition at line 61 of file ZElectronsSelector.cc.

Referenced by ZElectronsSelector().

vector<double> ZElectronsSelector::hOverECut

Definition at line 67 of file ZElectronsSelector.cc.

Referenced by operator()(), and ZElectronsSelector().

vector<int> ZElectronsSelector::missHits

Definition at line 63 of file ZElectronsSelector.cc.

Referenced by operator()(), and ZElectronsSelector().

vector<double> ZElectronsSelector::relCombIso

Definition at line 68 of file ZElectronsSelector.cc.

Referenced by operator()(), and ZElectronsSelector().

vector<double> ZElectronsSelector::sigmaIEtaIEtaCut

Definition at line 64 of file ZElectronsSelector.cc.

Referenced by operator()(), and ZElectronsSelector().

edm::EDGetTokenT<reco::GsfElectronCollection> ZElectronsSelector::theGsfEToken

Definition at line 54 of file ZElectronsSelector.cc.

edm::EDGetTokenT<double> ZElectronsSelector::theRhoToken

Definition at line 53 of file ZElectronsSelector.cc.

Referenced by newEvent().