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  {
17  contentTags_.emplace("rho",rhoTag);
18 
19  }
20 
21  result_type operator()(const reco::GsfElectronPtr&) const final;
22 
24  void getEventContent(const edm::EventBase&) final;
25 
26  double value(const reco::CandidatePtr& cand) const final;
27 
28  CandidateType candidateType() const final {
29  return ELECTRON;
30  }
31 
32 private:
35 };
36 
39  "GsfEleHadronicOverEMEnergyScaledCut");
40 
42  auto rho = cc.consumes<double>(contentTags_["rho"]);
43  contentTokens_.emplace("rho",rho);
44 }
45 
48 }
49 
50 
51 CutApplicatorBase::result_type GsfEleHadronicOverEMEnergyScaledCut::operator()(const reco::GsfElectronPtr& cand) const {
52 
53  const double rho = rhoHandle_.isValid() ? (float)(*rhoHandle_) : 0;
54  const float energy = cand->superCluster()->energy();
55  const float c0 = (std::abs(cand->superCluster()->position().eta()) < barrelCutOff_ ? barrelC0_ : endcapC0_);
56  const float cE = (std::abs(cand->superCluster()->position().eta()) < barrelCutOff_ ? barrelCE_ : endcapCE_);
57  const float cR = (std::abs(cand->superCluster()->position().eta()) < barrelCutOff_ ? barrelCr_ : endcapCr_);
58  return cand->hadronicOverEm() < c0 + cE/energy + cR*rho/energy;
59 }
60 
62  reco::GsfElectronPtr ele(cand);
63  return ele->hadronicOverEm();
64 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
result_type operator()(const reco::GsfElectronPtr &) const final
bool ev
std::unordered_map< std::string, edm::InputTag > contentTags_
float hadronicOverEm() const
Definition: GsfElectron.h:487
GsfEleHadronicOverEMEnergyScaledCut(const edm::ParameterSet &c)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isValid() const
Definition: HandleBase.h:74
void setConsumes(edm::ConsumesCollector &) final
double value(const reco::CandidatePtr &cand) const final
bool getByLabel(InputTag const &, Handle< T > &) const
Definition: EventBase.h:94
SuperClusterRef superCluster() const override
reference to a SuperCluster
Definition: GsfElectron.h:184
#define DEFINE_EDM_PLUGIN(factory, type, name)
void getEventContent(const edm::EventBase &) final