CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
GsfEleRelPFIsoScaledCut.cc
Go to the documentation of this file.
4 
6 public:
8 
9  result_type operator()(const reco::GsfElectronPtr&) const final;
10 
12  void getEventContent(const edm::EventBase&) final;
13 
14  double value(const reco::CandidatePtr& cand) const final;
15 
17 
18 private:
22 };
23 
25 
28  barrelC0_(c.getParameter<double>("barrelC0")),
29  endcapC0_(c.getParameter<double>("endcapC0")),
30  barrelCpt_(c.getParameter<double>("barrelCpt")),
31  endcapCpt_(c.getParameter<double>("endcapCpt")),
32  barrelCutOff_(c.getParameter<double>("barrelCutOff")),
33  effectiveAreas_((c.getParameter<edm::FileInPath>("effAreasConfigFile")).fullPath()) {
35  contentTags_.emplace("rho", rhoTag);
36 }
37 
39  auto rho = cc.consumes<double>(contentTags_["rho"]);
40  contentTokens_.emplace("rho", rho);
41 }
42 
45 }
46 
48  // Establish the cut value
49  double absEta = std::abs(cand->superCluster()->eta());
50 
51  const float C0 = (absEta < barrelCutOff_ ? barrelC0_ : endcapC0_);
52  const float Cpt = (absEta < barrelCutOff_ ? barrelCpt_ : endcapCpt_);
53  const float isoCut = C0 + Cpt / cand->pt();
54 
55  return value(cand) < isoCut;
56 }
57 
59  // Establish the cut value
60  reco::GsfElectronPtr ele(cand);
61  double absEta = std::abs(ele->superCluster()->eta());
62 
63  // Compute the combined isolation with effective area correction
64  auto pfIso = ele->pfIsolationVariables();
65  const float chad = pfIso.sumChargedHadronPt;
66  const float nhad = pfIso.sumNeutralHadronEt;
67  const float pho = pfIso.sumPhotonEt;
68  const float eA = effectiveAreas_.getEffectiveArea(absEta);
69  const float rho = rhoHandle_.isValid() ? (float)(*rhoHandle_) : 0; // std::max likes float arguments
70  const float iso = chad + std::max(0.0f, nhad + pho - rho * eA);
71  return iso / cand->pt();
72 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
const edm::EventSetup & c
std::unordered_map< std::string, edm::InputTag > contentTags_
bool ev
const float getEffectiveArea(float eta) const
GsfEleRelPFIsoScaledCut(const edm::ParameterSet &c)
CandidateType candidateType() const final
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::Handle< double > rhoHandle_
bool isValid() const
Definition: HandleBase.h:70
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
result_type operator()(const reco::GsfElectronPtr &) const final
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
void getEventContent(const edm::EventBase &) final
#define DEFINE_EDM_PLUGIN(factory, type, name)
double value(const reco::CandidatePtr &cand) const final