CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
GsfEleEffAreaPFIsoCut.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:
19  // Cut values
21  // Configuration
22  const float _ptCutOff;
23  const float _barrelCutOff;
25  // Effective area constants
27  // The rho
29 
30  constexpr static char rhoString_[] = "rho";
31 };
32 
33 constexpr char GsfEleEffAreaPFIsoCut::rhoString_[];
34 
36 
39  _isoCutEBLowPt(c.getParameter<double>("isoCutEBLowPt")),
40  _isoCutEBHighPt(c.getParameter<double>("isoCutEBHighPt")),
41  _isoCutEELowPt(c.getParameter<double>("isoCutEELowPt")),
42  _isoCutEEHighPt(c.getParameter<double>("isoCutEEHighPt")),
43  _ptCutOff(c.getParameter<double>("ptCutOff")),
44  _barrelCutOff(c.getParameter<double>("barrelCutOff")),
45  _isRelativeIso(c.getParameter<bool>("isRelativeIso")),
46  _effectiveAreas((c.getParameter<edm::FileInPath>("effAreasConfigFile")).fullPath()) {
48  contentTags_.emplace(rhoString_, rhoTag);
49 }
50 
52  auto rho = cc.consumes<double>(contentTags_[rhoString_]);
53  contentTokens_.emplace(rhoString_, rho);
54 }
55 
58 }
59 
61  // Establish the cut value
62  double absEta = std::abs(cand->superCluster()->eta());
63  const float isoCut = (cand->pt() < _ptCutOff ? (absEta < _barrelCutOff ? _isoCutEBLowPt : _isoCutEELowPt)
65 
66  // Compute the combined isolation with effective area correction
67  const reco::GsfElectron::PflowIsolationVariables& pfIso = cand->pfIsolationVariables();
68  const float chad = pfIso.sumChargedHadronPt;
69  const float nhad = pfIso.sumNeutralHadronEt;
70  const float pho = pfIso.sumPhotonEt;
71  const float eA = _effectiveAreas.getEffectiveArea(absEta);
72  const float rho = _rhoHandle.isValid() ? (float)(*_rhoHandle) : 0; // std::max likes float arguments
73  const float iso = chad + std::max(0.0f, nhad + pho - rho * eA);
74 
75  // Apply the cut and return the result
76  // Scale by pT if the relative isolation is requested but avoid division by 0
77  return iso < isoCut * (_isRelativeIso ? cand->pt() : 1.);
78 }
79 
81  reco::GsfElectronPtr ele(cand);
82  // Establish the cut value
83  double absEta = std::abs(ele->superCluster()->eta());
84 
85  // Compute the combined isolation with effective area correction
86  const reco::GsfElectron::PflowIsolationVariables& pfIso = ele->pfIsolationVariables();
87  const float chad = pfIso.sumChargedHadronPt;
88  const float nhad = pfIso.sumNeutralHadronEt;
89  const float pho = pfIso.sumPhotonEt;
90  float eA = _effectiveAreas.getEffectiveArea(absEta);
91  float rho = (float)(*_rhoHandle); // std::max likes float arguments
92  float iso = chad + std::max(0.0f, nhad + pho - rho * eA);
93 
94  // Divide by pT if the relative isolation is requested
95  if (_isRelativeIso)
96  iso /= ele->pt();
97 
98  // Apply the cut and return the result
99  return iso;
100 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
void getEventContent(const edm::EventBase &) final
const edm::EventSetup & c
std::unordered_map< std::string, edm::InputTag > contentTags_
bool ev
const float getEffectiveArea(float eta) const
float sumPhotonEt
sum pt of PF photons // old float photonIso ;
Definition: GsfElectron.h:665
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static constexpr char rhoString_[]
float sumNeutralHadronEt
sum pt of neutral hadrons // old float neutralHadronIso ;
Definition: GsfElectron.h:664
bool isValid() const
Definition: HandleBase.h:70
void setConsumes(edm::ConsumesCollector &) final
double value(const reco::CandidatePtr &cand) const final
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
bool getByLabel(InputTag const &, Handle< T > &) const
Definition: EventBase.h:92
result_type operator()(const reco::GsfElectronPtr &) const final
#define DEFINE_EDM_PLUGIN(factory, type, name)
GsfEleEffAreaPFIsoCut(const edm::ParameterSet &c)
edm::Handle< double > _rhoHandle
CandidateType candidateType() const final
float sumChargedHadronPt
sum-pt of charged Hadron // old float chargedHadronIso ;
Definition: GsfElectron.h:663