CMS 3D CMS Logo

GsfEleRelPFIsoScaledCut.cc
Go to the documentation of this file.
4 
5 
7  public:
9 
10  result_type operator()(const reco::GsfElectronPtr&) const final;
11 
13  void getEventContent(const edm::EventBase&) final;
14 
15  double value(const reco::CandidatePtr& cand) const final;
16 
17  CandidateType candidateType() const final {
18  return ELECTRON;
19  }
20 
21  private:
25 };
26 
28 
31  barrelC0_(c.getParameter<double>("barrelC0")),
32  endcapC0_(c.getParameter<double>("endcapC0")),
33  barrelCpt_(c.getParameter<double>("barrelCpt")),
34  endcapCpt_(c.getParameter<double>("endcapCpt")),
35  barrelCutOff_(c.getParameter<double>("barrelCutOff")),
36  effectiveAreas_((c.getParameter<edm::FileInPath>("effAreasConfigFile")).fullPath())
37 {
39  contentTags_.emplace("rho",rhoTag);
40 }
41 
43  auto rho = cc.consumes<double>(contentTags_["rho"]);
44  contentTokens_.emplace("rho", rho);
45 }
46 
49 }
50 
51 CutApplicatorBase::result_type GsfEleRelPFIsoScaledCut::operator()(const reco::GsfElectronPtr& cand) const {
52  // Establish the cut value
53  double absEta = std::abs(cand->superCluster()->eta());
54 
55  const float C0 = (absEta < barrelCutOff_ ? barrelC0_ : endcapC0_);
56  const float Cpt = (absEta < barrelCutOff_ ? barrelCpt_ : endcapCpt_);
57  const float isoCut = C0+Cpt/cand->pt();
58 
59  return value(cand) < isoCut;
60 }
61 
63 
64  // Establish the cut value
65  reco::GsfElectronPtr ele(cand);
66  double absEta = std::abs(ele->superCluster()->eta());
67 
68  // Compute the combined isolation with effective area correction
69  auto pfIso = ele->pfIsolationVariables();
70  const float chad = pfIso.sumChargedHadronPt;
71  const float nhad = pfIso.sumNeutralHadronEt;
72  const float pho = pfIso.sumPhotonEt;
73  const float eA = effectiveAreas_.getEffectiveArea(absEta);
74  const float rho = rhoHandle_.isValid() ? (float)(*rhoHandle_) : 0; // std::max likes float arguments
75  const float iso = chad + std::max(0.0f, nhad + pho - rho*eA);
76  return iso/cand->pt();
77 }
const PflowIsolationVariables & pfIsolationVariables() const
Definition: GsfElectron.h:677
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
double pt() const final
transverse momentum
bool ev
const float getEffectiveArea(float eta) const
std::unordered_map< std::string, edm::InputTag > contentTags_
GsfEleRelPFIsoScaledCut(const edm::ParameterSet &c)
result_type operator()(const reco::GsfElectronPtr &) const final
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
edm::Handle< double > rhoHandle_
bool isValid() const
Definition: HandleBase.h:74
void setConsumes(edm::ConsumesCollector &) final
CandidateType candidateType() const final
double value(const reco::CandidatePtr &cand) const final
bool getByLabel(InputTag const &, Handle< T > &) const
Definition: EventBase.h:94
HLT enums.
void getEventContent(const edm::EventBase &) final
SuperClusterRef superCluster() const override
reference to a SuperCluster
Definition: GsfElectron.h:184
#define DEFINE_EDM_PLUGIN(factory, type, name)
float sumChargedHadronPt
sum-pt of charged Hadron // old float chargedHadronIso ;
Definition: GsfElectron.h:629