CMS 3D CMS Logo

GsfEleDeltaBetaIsoCutStandalone.cc
Go to the documentation of this file.
6 
7 
9 public:
11 
12  result_type operator()(const reco::GsfElectronPtr&) const final;
13 
14  double value(const reco::CandidatePtr& cand) const final;
15 
16  CandidateType candidateType() const final {
17  return ELECTRON;
18  }
19 
20 private:
23  const bool _relativeIso;
24 
25 };
26 
29  "GsfEleDeltaBetaIsoCutStandalone");
30 
33  _isoCutEBLowPt(c.getParameter<double>("isoCutEBLowPt")),
34  _isoCutEBHighPt(c.getParameter<double>("isoCutEBHighPt")),
35  _isoCutEELowPt(c.getParameter<double>("isoCutEELowPt")),
36  _isoCutEEHighPt(c.getParameter<double>("isoCutEEHighPt")),
37  _deltaBetaConstant(c.getParameter<double>("deltaBetaConstant")),
38  _ptCutOff(c.getParameter<double>("ptCutOff")),
39  _barrelCutOff(c.getParameter<double>("barrelCutOff")),
40  _relativeIso(c.getParameter<bool>("isRelativeIso")) {
41 }
42 
43 CutApplicatorBase::result_type
46  const float isoCut =
47  ( cand->p4().pt() < _ptCutOff ?
48  ( std::abs(cand->superCluster()->position().eta()) < _barrelCutOff ?
50  ( std::abs(cand->superCluster()->position().eta()) < _barrelCutOff ?
53  cand->pfIsolationVariables();
54  const float chad = pfIso.sumChargedHadronPt;
55  const float nhad = pfIso.sumNeutralHadronEt;
56  const float pho = pfIso.sumPhotonEt;
57  const float puchad = pfIso.sumPUPt;
58  float iso = chad + std::max(0.0f, nhad + pho - _deltaBetaConstant*puchad);
59  if( _relativeIso ) iso /= cand->p4().pt();
60  return iso < isoCut;
61 }
62 
64 value(const reco::CandidatePtr& cand) const {
65  reco::GsfElectronPtr ele(cand);
67  ele->pfIsolationVariables();
68  const float chad = pfIso.sumChargedHadronPt;
69  const float nhad = pfIso.sumNeutralHadronEt;
70  const float pho = pfIso.sumPhotonEt;
71  const float puchad = pfIso.sumPUPt;
72  float iso = chad + std::max(0.0f, nhad + pho - _deltaBetaConstant*puchad);
73  if( _relativeIso ) iso /= cand->p4().pt();
74  return iso;
75 }
const PflowIsolationVariables & pfIsolationVariables() const
Definition: GsfElectron.h:670
const LorentzVector & p4(P4Kind kind) const
Definition: GsfElectron.cc:225
float sumPUPt
sum pt of charged Particles not from PV (for Pu corrections)
Definition: GsfElectron.h:632
GsfEleDeltaBetaIsoCutStandalone(const edm::ParameterSet &c)
float sumPhotonEt
sum pt of PF photons // old float photonIso ;
Definition: GsfElectron.h:627
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
float sumNeutralHadronEt
sum pt of neutral hadrons // old float neutralHadronIso ;
Definition: GsfElectron.h:626
double value(const reco::CandidatePtr &cand) const final
result_type operator()(const reco::GsfElectronPtr &) const 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:625