CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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")) {
16  contentTags_.emplace("rho", rhoTag);
17  }
18 
20 
21  void setConsumes(edm::ConsumesCollector&) final;
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 
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 
54  reco::GsfElectronPtr ele(cand);
55  return ele->hadronicOverEm();
56 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
const edm::EventSetup & c
result_type operator()(const reco::GsfElectronPtr &) const final
std::unordered_map< std::string, edm::InputTag > contentTags_
double value(const reco::CandidatePtr &cand) const final
bool ev
GsfEleHadronicOverEMEnergyScaledCut(const edm::ParameterSet &c)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isValid() const
Definition: HandleBase.h:70
edm::Ptr< Candidate > CandidatePtr
persistent reference to an object in a collection of Candidate objects
Definition: CandidateFwd.h:25
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
void setConsumes(edm::ConsumesCollector &) final
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
bool getByLabel(InputTag const &, Handle< T > &) const
Definition: EventBase.h:92
#define DEFINE_EDM_PLUGIN(factory, type, name)
void getEventContent(const edm::EventBase &) final