CMS 3D CMS Logo

GsfEleHadronicOverEMEnergyScaledCut.cc
Go to the documentation of this file.
3 
5 public:
8  barrelC0_(c.getParameter<double>("barrelC0")),
9  barrelCE_(c.getParameter<double>("barrelCE")),
10  barrelCr_(c.getParameter<double>("barrelCr")),
11  endcapC0_(c.getParameter<double>("endcapC0")),
12  endcapCE_(c.getParameter<double>("endcapCE")),
13  endcapCr_(c.getParameter<double>("endcapCr")),
14  barrelCutOff_(c.getParameter<double>("barrelCutOff")) {
15  edm::InputTag rhoTag = c.getParameter<edm::InputTag>("rho");
16  contentTags_.emplace("rho", rhoTag);
17  }
18 
19  result_type operator()(const reco::GsfElectronPtr&) const final;
20 
22  void getEventContent(const edm::EventBase&) final;
23 
24  double value(const reco::CandidatePtr& cand) const final;
25 
27 
28 private:
31 };
32 
33 DEFINE_EDM_PLUGIN(CutApplicatorFactory, GsfEleHadronicOverEMEnergyScaledCut, "GsfEleHadronicOverEMEnergyScaledCut");
34 
36  auto rho = cc.consumes<double>(contentTags_["rho"]);
37  contentTokens_.emplace("rho", rho);
38 }
39 
41  ev.getByLabel(contentTags_["rho"], rhoHandle_);
42 }
43 
45  const double rho = rhoHandle_.isValid() ? (float)(*rhoHandle_) : 0;
46  const float energy = cand->superCluster()->energy();
47  const float c0 = (std::abs(cand->superCluster()->position().eta()) < barrelCutOff_ ? barrelC0_ : endcapC0_);
48  const float cE = (std::abs(cand->superCluster()->position().eta()) < barrelCutOff_ ? barrelCE_ : endcapCE_);
49  const float cR = (std::abs(cand->superCluster()->position().eta()) < barrelCutOff_ ? barrelCr_ : endcapCr_);
50  return cand->hadronicOverEm() < c0 + cE / energy + cR * rho / energy;
51 }
52 
55  return ele->hadronicOverEm();
56 }
GsfEleHadronicOverEMEnergyScaledCut::endcapCr_
const float endcapCr_
Definition: GsfEleHadronicOverEMEnergyScaledCut.cc:29
GsfEleHadronicOverEMEnergyScaledCut::operator()
result_type operator()(const reco::GsfElectronPtr &) const final
Definition: GsfEleHadronicOverEMEnergyScaledCut.cc:44
dqmMemoryStats.float
float
Definition: dqmMemoryStats.py:127
GsfEleHadronicOverEMEnergyScaledCut::barrelCr_
const float barrelCr_
Definition: GsfEleHadronicOverEMEnergyScaledCut.cc:29
GsfEleHadronicOverEMEnergyScaledCut::barrelCE_
const float barrelCE_
Definition: GsfEleHadronicOverEMEnergyScaledCut.cc:29
GsfEleHadronicOverEMEnergyScaledCut::value
double value(const reco::CandidatePtr &cand) const final
Definition: GsfEleHadronicOverEMEnergyScaledCut.cc:53
reco::GsfElectron::hadronicOverEm
float hadronicOverEm() const
Definition: GsfElectron.h:475
watchdog.const
const
Definition: watchdog.py:83
edm::Handle< double >
CutApplicatorWithEventContentBase::contentTags_
std::unordered_map< std::string, edm::InputTag > contentTags_
Definition: CutApplicatorWithEventContentBase.h:35
HLT_FULL_cff.rhoTag
rhoTag
Definition: HLT_FULL_cff.py:14988
candidate_functions::CandidateCut::result_type
bool result_type
Definition: CandidateCut.h:11
GsfEleHadronicOverEMEnergyScaledCut::barrelCutOff_
const float barrelCutOff_
Definition: GsfEleHadronicOverEMEnergyScaledCut.cc:29
CutApplicatorWithEventContentBase.h
GsfEleHadronicOverEMEnergyScaledCut
Definition: GsfEleHadronicOverEMEnergyScaledCut.cc:4
HCALHighEnergyHPDFilter_cfi.energy
energy
Definition: HCALHighEnergyHPDFilter_cfi.py:5
GsfEleHadronicOverEMEnergyScaledCut::barrelC0_
const float barrelC0_
Definition: GsfEleHadronicOverEMEnergyScaledCut.cc:29
GsfElectron.h
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition: PluginFactory.h:124
DDAxes::rho
edm::ParameterSet
Definition: ParameterSet.h:47
GsfEleHadronicOverEMEnergyScaledCut::setConsumes
void setConsumes(edm::ConsumesCollector &) final
Definition: GsfEleHadronicOverEMEnergyScaledCut.cc:35
CutApplicatorBase::CandidateType
CandidateType
Definition: CutApplicatorBase.h:47
edmplugin::PluginFactory
Definition: PluginFactory.h:34
GsfEleHadronicOverEMEnergyScaledCut::GsfEleHadronicOverEMEnergyScaledCut
GsfEleHadronicOverEMEnergyScaledCut(const edm::ParameterSet &c)
Definition: GsfEleHadronicOverEMEnergyScaledCut.cc:6
cand
Definition: decayParser.h:32
GsfEleHadronicOverEMEnergyScaledCut::candidateType
CandidateType candidateType() const final
Definition: GsfEleHadronicOverEMEnergyScaledCut.cc:26
GsfEleHadronicOverEMEnergyScaledCut::endcapC0_
const float endcapC0_
Definition: GsfEleHadronicOverEMEnergyScaledCut.cc:29
cc
edm::Ptr< reco::GsfElectron >
ev
bool ev
Definition: Hydjet2Hadronizer.cc:95
CutApplicatorWithEventContentBase::contentTokens_
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
Definition: CutApplicatorWithEventContentBase.h:40
edm::EventBase
Definition: EventBase.h:46
fftjetpileupestimator_calo_uncalib_cfi.c0
c0
Definition: fftjetpileupestimator_calo_uncalib_cfi.py:8
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
GsfEleHadronicOverEMEnergyScaledCut::getEventContent
void getEventContent(const edm::EventBase &) final
Definition: GsfEleHadronicOverEMEnergyScaledCut.cc:40
c
auto & c
Definition: CAHitNtupletGeneratorKernelsImpl.h:46
edm::HandleBase::isValid
bool isValid() const
Definition: HandleBase.h:70
CutApplicatorBase::ELECTRON
Definition: CutApplicatorBase.h:47
edm::InputTag
Definition: InputTag.h:15
edm::ConsumesCollector
Definition: ConsumesCollector.h:45
CutApplicatorWithEventContentBase
Definition: CutApplicatorWithEventContentBase.h:19
GsfEleHadronicOverEMEnergyScaledCut::endcapCE_
const float endcapCE_
Definition: GsfEleHadronicOverEMEnergyScaledCut.cc:29
GsfEleHadronicOverEMEnergyScaledCut::rhoHandle_
edm::Handle< double > rhoHandle_
Definition: GsfEleHadronicOverEMEnergyScaledCut.cc:30