CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GsfEleDeltaBetaIsoCut.cc
Go to the documentation of this file.
6 
7 
9 public:
11 
12  result_type operator()(const reco::GsfElectronPtr&) const override final;
13 
14  void setConsumes(edm::ConsumesCollector&) override final;
15  void getEventContent(const edm::EventBase&) override final;
16 
17  CandidateType candidateType() const override final {
18  return ELECTRON;
19  }
20 
21 private:
26 };
27 
30  "GsfEleDeltaBetaIsoCut");
31 
34  _isoCutEBLowPt(c.getParameter<double>("isoCutEBLowPt")),
35  _isoCutEBHighPt(c.getParameter<double>("isoCutEBHighPt")),
36  _isoCutEELowPt(c.getParameter<double>("isoCutEELowPt")),
37  _isoCutEEHighPt(c.getParameter<double>("isoCutEEHighPt")),
38  _deltaBetaConstant(c.getParameter<double>("deltaBetaConstant")),
39  _ptCutOff(c.getParameter<double>("ptCutOff")),
40  _barrelCutOff(c.getParameter<double>("barrelCutOff")),
41  _relativeIso(c.getParameter<bool>("isRelativeIso")) {
42  edm::InputTag chadtag = c.getParameter<edm::InputTag>("chargedHadronIsoDR03Src");
43  edm::InputTag nhadtag = c.getParameter<edm::InputTag>("neutralHadronIsoDR03Src");
44  edm::InputTag phtag = c.getParameter<edm::InputTag>("photonIsoDR03Src");
45  edm::InputTag PUchadtag = c.getParameter<edm::InputTag>("puChargedHadronIsoDR03Src");
46  contentTags_.emplace("chad",chadtag);
47  contentTags_.emplace("nhad",nhadtag);
48  contentTags_.emplace("pho",phtag);
49  contentTags_.emplace("puchad",PUchadtag);
50 }
51 
53  auto chad =
55  contentTokens_.emplace("chad",chad);
56  auto nhad =
58  contentTokens_.emplace("nhad",nhad);
59  auto pho =
61  contentTokens_.emplace("pho",pho);
62  auto puchad =
64  contentTokens_.emplace("puchad",puchad);
65 }
66 
68  ev.getByLabel(contentTags_["chad"],_chad_iso);
69  ev.getByLabel(contentTags_["nhad"],_nhad_iso);
70  ev.getByLabel(contentTags_["pho"],_ph_iso);
71  ev.getByLabel(contentTags_["puchad"],_PUchad_iso);
72 }
73 
74 CutApplicatorBase::result_type
76 operator()(const reco::GsfElectronPtr& cand) const{
77  const float isoCut =
78  ( cand->p4().pt() < _ptCutOff ?
79  ( std::abs(cand->superCluster()->position().eta()) < _barrelCutOff ?
81  ( std::abs(cand->superCluster()->position().eta()) < _barrelCutOff ?
83  const float chad = (*_chad_iso)[cand];
84  const float nhad = (*_nhad_iso)[cand];
85  const float pho = (*_ph_iso)[cand];
86  const float puchad = (*_PUchad_iso)[cand];
87  float iso = chad + std::max(0.0f, nhad + pho - _deltaBetaConstant*puchad);
88  if( _relativeIso ) iso /= cand->p4().pt();
89  return iso < isoCut;
90 }
edm::Handle< edm::ValueMap< float > > _nhad_iso
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
void getEventContent(const edm::EventBase &) overridefinal
void setConsumes(edm::ConsumesCollector &) overridefinal
edm::Handle< edm::ValueMap< float > > _chad_iso
edm::Handle< edm::ValueMap< float > > _ph_iso
result_type operator()(const reco::GsfElectronPtr &) const overridefinal
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
CandidateType candidateType() const overridefinal
std::unordered_map< std::string, edm::InputTag > contentTags_
const T & max(const T &a, const T &b)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
edm::Handle< edm::ValueMap< float > > _PUchad_iso
string const
Definition: compareJSON.py:14
bool getByLabel(InputTag const &, Handle< T > &) const
Definition: EventBase.h:87
#define DEFINE_EDM_PLUGIN(factory, type, name)
GsfEleDeltaBetaIsoCut(const edm::ParameterSet &c)