CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes | Static Private Attributes
GsfEleEffAreaPFIsoCut Class Reference
Inheritance diagram for GsfEleEffAreaPFIsoCut:
CutApplicatorWithEventContentBase CutApplicatorBase candidate_functions::CandidateCut

Public Member Functions

CandidateType candidateType () const final
 
void getEventContent (const edm::EventBase &) final
 
 GsfEleEffAreaPFIsoCut (const edm::ParameterSet &c)
 
result_type operator() (const reco::GsfElectronPtr &) const final
 
void setConsumes (edm::ConsumesCollector &) final
 
double value (const reco::CandidatePtr &cand) const final
 
- Public Member Functions inherited from CutApplicatorWithEventContentBase
 CutApplicatorWithEventContentBase ()
 
 CutApplicatorWithEventContentBase (const edm::ParameterSet &c)
 
 CutApplicatorWithEventContentBase (const CutApplicatorWithEventContentBase &)=delete
 
CutApplicatorWithEventContentBaseoperator= (const CutApplicatorWithEventContentBase &)=delete
 
 ~CutApplicatorWithEventContentBase () override
 Destructor. More...
 
- Public Member Functions inherited from CutApplicatorBase
virtual result_type asCandidate (const argument_type &) const
 
 CutApplicatorBase ()
 
 CutApplicatorBase (const edm::ParameterSet &c)
 
 CutApplicatorBase (const CutApplicatorBase &)=delete
 
const std::string & name () const override
 
result_type operator() (const argument_type &) const final
 
virtual result_type operator() (const pat::ElectronPtr &) const
 
virtual result_type operator() (const reco::PhotonPtr &) const
 
virtual result_type operator() (const pat::PhotonPtr &) const
 
virtual result_type operator() (const reco::MuonPtr &) const
 
virtual result_type operator() (const pat::MuonPtr &) const
 
virtual result_type operator() (const reco::PFTauPtr &) const
 
virtual result_type operator() (const pat::TauPtr &) const
 
CutApplicatorBaseoperator= (const CutApplicatorBase &)=delete
 
 ~CutApplicatorBase () override
 Destructor. More...
 
- Public Member Functions inherited from candidate_functions::CandidateCut
 CandidateCut ()
 
virtual ~CandidateCut ()
 

Private Attributes

const float _barrelCutOff
 
EffectiveAreas _effectiveAreas
 
const float _isoCutEBHighPt
 
const float _isoCutEBLowPt
 
const float _isoCutEEHighPt
 
const float _isoCutEELowPt
 
bool _isRelativeIso
 
const float _ptCutOff
 
edm::Handle< double > _rhoHandle
 

Static Private Attributes

static constexpr char rhoString_ [] = "rho"
 

Additional Inherited Members

- Public Types inherited from CutApplicatorBase
enum  CandidateType {
  NONE, ELECTRON, MUON, PHOTON,
  TAU, PATELECTRON, PATMUON, PATPHOTON,
  PATTAU
}
 
- Public Types inherited from candidate_functions::CandidateCut
using argument_type = reco::CandidatePtr
 
using result_type = bool
 
- Protected Attributes inherited from CutApplicatorWithEventContentBase
std::unordered_map< std::string, edm::InputTagcontentTags_
 
std::unordered_map< std::string, edm::EDGetTokencontentTokens_
 

Detailed Description

Definition at line 5 of file GsfEleEffAreaPFIsoCut.cc.

Constructor & Destructor Documentation

◆ GsfEleEffAreaPFIsoCut()

GsfEleEffAreaPFIsoCut::GsfEleEffAreaPFIsoCut ( const edm::ParameterSet c)

Definition at line 37 of file GsfEleEffAreaPFIsoCut.cc.

References HltBtagPostValidation_cff::c, CutApplicatorWithEventContentBase::contentTags_, rhoString_, and HLT_2024v14_cff::rhoTag.

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()) {
47  edm::InputTag rhoTag = c.getParameter<edm::InputTag>("rho");
48  contentTags_.emplace(rhoString_, rhoTag);
49 }
std::unordered_map< std::string, edm::InputTag > contentTags_
static constexpr char rhoString_[]

Member Function Documentation

◆ candidateType()

CandidateType GsfEleEffAreaPFIsoCut::candidateType ( ) const
inlinefinalvirtual

Reimplemented from CutApplicatorBase.

Definition at line 16 of file GsfEleEffAreaPFIsoCut.cc.

References CutApplicatorBase::ELECTRON.

◆ getEventContent()

void GsfEleEffAreaPFIsoCut::getEventContent ( const edm::EventBase ev)
finalvirtual

◆ operator()()

CutApplicatorBase::result_type GsfEleEffAreaPFIsoCut::operator() ( const reco::GsfElectronPtr cand) const
finalvirtual

Reimplemented from CutApplicatorBase.

Definition at line 60 of file GsfEleEffAreaPFIsoCut.cc.

References _barrelCutOff, _effectiveAreas, _isoCutEBHighPt, _isoCutEBLowPt, _isoCutEEHighPt, _isoCutEELowPt, _isRelativeIso, _ptCutOff, _rhoHandle, funct::abs(), f, nano_mu_digi_cff::float, EffectiveAreas::getEffectiveArea(), genparticles_cff::iso, edm::HandleBase::isValid(), SiStripPI::max, rho, reco::GsfElectron::PflowIsolationVariables::sumChargedHadronPt, reco::GsfElectron::PflowIsolationVariables::sumNeutralHadronEt, and reco::GsfElectron::PflowIsolationVariables::sumPhotonEt.

60  {
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 }
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
double f[11][100]
float sumNeutralHadronEt
sum pt of neutral hadrons // old float neutralHadronIso ;
Definition: GsfElectron.h:664
const float getEffectiveArea(float eta) const
bool isValid() const
Definition: HandleBase.h:70
edm::Handle< double > _rhoHandle
float sumChargedHadronPt
sum-pt of charged Hadron // old float chargedHadronIso ;
Definition: GsfElectron.h:663

◆ setConsumes()

void GsfEleEffAreaPFIsoCut::setConsumes ( edm::ConsumesCollector cc)
finalvirtual

Implements CutApplicatorWithEventContentBase.

Definition at line 51 of file GsfEleEffAreaPFIsoCut.cc.

References gpuPixelDoublets::cc, CutApplicatorWithEventContentBase::contentTags_, CutApplicatorWithEventContentBase::contentTokens_, rho, and rhoString_.

51  {
52  auto rho = cc.consumes<double>(contentTags_[rhoString_]);
53  contentTokens_.emplace(rhoString_, rho);
54 }
std::unordered_map< std::string, edm::InputTag > contentTags_
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
static constexpr char rhoString_[]
std::unordered_map< std::string, edm::EDGetToken > contentTokens_

◆ value()

double GsfEleEffAreaPFIsoCut::value ( const reco::CandidatePtr cand) const
finalvirtual

Implements candidate_functions::CandidateCut.

Definition at line 80 of file GsfEleEffAreaPFIsoCut.cc.

References _effectiveAreas, _isRelativeIso, _rhoHandle, funct::abs(), f, nano_mu_digi_cff::float, EffectiveAreas::getEffectiveArea(), genparticles_cff::iso, SiStripPI::max, reco::GsfElectron::pfIsolationVariables(), reco::LeafCandidate::pt(), rho, reco::GsfElectron::PflowIsolationVariables::sumChargedHadronPt, and reco::GsfElectron::superCluster().

80  {
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 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
const float getEffectiveArea(float eta) const
edm::Handle< double > _rhoHandle
float sumChargedHadronPt
sum-pt of charged Hadron // old float chargedHadronIso ;
Definition: GsfElectron.h:663

Member Data Documentation

◆ _barrelCutOff

const float GsfEleEffAreaPFIsoCut::_barrelCutOff
private

Definition at line 23 of file GsfEleEffAreaPFIsoCut.cc.

Referenced by operator()().

◆ _effectiveAreas

EffectiveAreas GsfEleEffAreaPFIsoCut::_effectiveAreas
private

Definition at line 26 of file GsfEleEffAreaPFIsoCut.cc.

Referenced by operator()(), and value().

◆ _isoCutEBHighPt

const float GsfEleEffAreaPFIsoCut::_isoCutEBHighPt
private

Definition at line 20 of file GsfEleEffAreaPFIsoCut.cc.

Referenced by operator()().

◆ _isoCutEBLowPt

const float GsfEleEffAreaPFIsoCut::_isoCutEBLowPt
private

Definition at line 20 of file GsfEleEffAreaPFIsoCut.cc.

Referenced by operator()().

◆ _isoCutEEHighPt

const float GsfEleEffAreaPFIsoCut::_isoCutEEHighPt
private

Definition at line 20 of file GsfEleEffAreaPFIsoCut.cc.

Referenced by operator()().

◆ _isoCutEELowPt

const float GsfEleEffAreaPFIsoCut::_isoCutEELowPt
private

Definition at line 20 of file GsfEleEffAreaPFIsoCut.cc.

Referenced by operator()().

◆ _isRelativeIso

bool GsfEleEffAreaPFIsoCut::_isRelativeIso
private

Definition at line 24 of file GsfEleEffAreaPFIsoCut.cc.

Referenced by operator()(), and value().

◆ _ptCutOff

const float GsfEleEffAreaPFIsoCut::_ptCutOff
private

Definition at line 22 of file GsfEleEffAreaPFIsoCut.cc.

Referenced by operator()().

◆ _rhoHandle

edm::Handle<double> GsfEleEffAreaPFIsoCut::_rhoHandle
private

Definition at line 28 of file GsfEleEffAreaPFIsoCut.cc.

Referenced by getEventContent(), operator()(), and value().

◆ rhoString_

constexpr char GsfEleEffAreaPFIsoCut::rhoString_ = "rho"
staticprivate

Definition at line 30 of file GsfEleEffAreaPFIsoCut.cc.

Referenced by getEventContent(), GsfEleEffAreaPFIsoCut(), and setConsumes().