CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 overridefinal
 
void getEventContent (const edm::EventBase &) overridefinal
 
 GsfEleEffAreaPFIsoCut (const edm::ParameterSet &c)
 
result_type operator() (const reco::GsfElectronPtr &) const overridefinal
 
void setConsumes (edm::ConsumesCollector &) overridefinal
 
double value (const reco::CandidatePtr &cand) const overridefinal
 
- Public Member Functions inherited from CutApplicatorWithEventContentBase
 CutApplicatorWithEventContentBase (const edm::ParameterSet &c)
 
 CutApplicatorWithEventContentBase (const CutApplicatorWithEventContentBase &)=delete
 
CutApplicatorWithEventContentBaseoperator= (const CutApplicatorWithEventContentBase &)=delete
 
virtual ~CutApplicatorWithEventContentBase ()
 Destructor. More...
 
- Public Member Functions inherited from CutApplicatorBase
virtual result_type asCandidate (const argument_type &) const
 
 CutApplicatorBase (const edm::ParameterSet &c)
 
 CutApplicatorBase (const CutApplicatorBase &)=delete
 
const std::string & name () const
 
virtual 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
 
virtual ~CutApplicatorBase ()
 Destructor. More...
 
- Public Member Functions inherited from candidate_functions::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 char rhoString_ [] = "rho"
 

Additional Inherited Members

- Public Types inherited from CutApplicatorBase
enum  CandidateType {
  NONE, ELECTRON, MUON, PHOTON,
  TAU, PATELECTRON, PATMUON, PATPHOTON,
  PATTAU
}
 
- Protected Attributes inherited from CutApplicatorWithEventContentBase
std::unordered_map
< std::string, edm::InputTag
contentTags_
 
std::unordered_map
< std::string, edm::EDGetToken
contentTokens_
 

Detailed Description

Definition at line 6 of file GsfEleEffAreaPFIsoCut.cc.

Constructor & Destructor Documentation

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

Definition at line 42 of file GsfEleEffAreaPFIsoCut.cc.

References CutApplicatorWithEventContentBase::contentTags_, edm::ParameterSet::getParameter(), and rhoString_.

42  :
44  _isoCutEBLowPt(c.getParameter<double>("isoCutEBLowPt")),
45  _isoCutEBHighPt(c.getParameter<double>("isoCutEBHighPt")),
46  _isoCutEELowPt(c.getParameter<double>("isoCutEELowPt")),
47  _isoCutEEHighPt(c.getParameter<double>("isoCutEEHighPt")),
48  _ptCutOff(c.getParameter<double>("ptCutOff")),
49  _barrelCutOff(c.getParameter<double>("barrelCutOff")),
50  _isRelativeIso(c.getParameter<bool>("isRelativeIso")),
51  _effectiveAreas( (c.getParameter<edm::FileInPath>("effAreasConfigFile")).fullPath())
52 {
53 
54  edm::InputTag rhoTag = c.getParameter<edm::InputTag>("rho");
55  contentTags_.emplace(rhoString_,rhoTag);
56 
57 }
T getParameter(std::string const &) const
CutApplicatorWithEventContentBase(const edm::ParameterSet &c)
std::unordered_map< std::string, edm::InputTag > contentTags_

Member Function Documentation

CandidateType GsfEleEffAreaPFIsoCut::candidateType ( ) const
inlinefinaloverridevirtual

Reimplemented from CutApplicatorBase.

Definition at line 17 of file GsfEleEffAreaPFIsoCut.cc.

References CutApplicatorBase::ELECTRON.

17  {
18  return ELECTRON;
19  }
void GsfEleEffAreaPFIsoCut::getEventContent ( const edm::EventBase ev)
finaloverridevirtual

Implements CutApplicatorWithEventContentBase.

Definition at line 65 of file GsfEleEffAreaPFIsoCut.cc.

References _rhoHandle, CutApplicatorWithEventContentBase::contentTags_, edm::EventBase::getByLabel(), and rhoString_.

65  {
66 
68 }
std::unordered_map< std::string, edm::InputTag > contentTags_
bool getByLabel(InputTag const &, Handle< T > &) const
Definition: EventBase.h:92
edm::Handle< double > _rhoHandle
CutApplicatorBase::result_type GsfEleEffAreaPFIsoCut::operator() ( const reco::GsfElectronPtr cand) const
finaloverridevirtual

Reimplemented from CutApplicatorBase.

Definition at line 72 of file GsfEleEffAreaPFIsoCut.cc.

References _barrelCutOff, _effectiveAreas, _isoCutEBHighPt, _isoCutEBLowPt, _isoCutEEHighPt, _isoCutEELowPt, _isRelativeIso, _ptCutOff, _rhoHandle, funct::abs(), f, EffectiveAreas::getEffectiveArea(), bookConverter::max, rho, reco::GsfElectron::PflowIsolationVariables::sumChargedHadronPt, reco::GsfElectron::PflowIsolationVariables::sumNeutralHadronEt, and reco::GsfElectron::PflowIsolationVariables::sumPhotonEt.

72  {
73 
74  // Establish the cut value
75  double absEta = std::abs(cand->superCluster()->eta());
76  const float isoCut =
77  ( cand->pt() < _ptCutOff ?
79  :
81 
82  // Compute the combined isolation with effective area correction
84  cand->pfIsolationVariables();
85  const float chad = pfIso.sumChargedHadronPt;
86  const float nhad = pfIso.sumNeutralHadronEt;
87  const float pho = pfIso.sumPhotonEt;
88  float eA = _effectiveAreas.getEffectiveArea( absEta );
89  float rho = (float)(*_rhoHandle); // std::max likes float arguments
90  float iso = chad + std::max(0.0f, nhad + pho - rho*eA);
91 
92  // Divide by pT if the relative isolation is requested
93  if( _isRelativeIso )
94  iso /= cand->pt();
95 
96  // Apply the cut and return the result
97  return iso < isoCut;
98 }
Definition: DDAxes.h:10
const float getEffectiveArea(float eta) const
float sumPhotonEt
sum pt of PF photons // old float photonIso ;
Definition: GsfElectron.h:561
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:560
edm::Handle< double > _rhoHandle
float sumChargedHadronPt
sum-pt of charged Hadron // old float chargedHadronIso ;
Definition: GsfElectron.h:559
void GsfEleEffAreaPFIsoCut::setConsumes ( edm::ConsumesCollector cc)
finaloverridevirtual

Implements CutApplicatorWithEventContentBase.

Definition at line 59 of file GsfEleEffAreaPFIsoCut.cc.

References edm::ConsumesCollector::consumes(), CutApplicatorWithEventContentBase::contentTags_, CutApplicatorWithEventContentBase::contentTokens_, rho, and rhoString_.

59  {
60 
61  auto rho = cc.consumes<double>(contentTags_[rhoString_]);
62  contentTokens_.emplace(rhoString_, rho);
63 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
Definition: DDAxes.h:10
std::unordered_map< std::string, edm::InputTag > contentTags_
double GsfEleEffAreaPFIsoCut::value ( const reco::CandidatePtr cand) const
finaloverridevirtual

Implements candidate_functions::CandidateCut.

Definition at line 100 of file GsfEleEffAreaPFIsoCut.cc.

References _effectiveAreas, _isRelativeIso, _rhoHandle, funct::abs(), f, EffectiveAreas::getEffectiveArea(), bookConverter::max, and rho.

100  {
101  reco::GsfElectronPtr ele(cand);
102  // Establish the cut value
103  double absEta = std::abs(ele->superCluster()->eta());
104 
105  // Compute the combined isolation with effective area correction
107  ele->pfIsolationVariables();
108  const float chad = pfIso.sumChargedHadronPt;
109  const float nhad = pfIso.sumNeutralHadronEt;
110  const float pho = pfIso.sumPhotonEt;
111  float eA = _effectiveAreas.getEffectiveArea( absEta );
112  float rho = (float)(*_rhoHandle); // std::max likes float arguments
113  float iso = chad + std::max(0.0f, nhad + pho - rho*eA);
114 
115  // Divide by pT if the relative isolation is requested
116  if( _isRelativeIso )
117  iso /= ele->pt();
118 
119  // Apply the cut and return the result
120  return iso;
121 }
Definition: DDAxes.h:10
const float getEffectiveArea(float eta) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
edm::Handle< double > _rhoHandle
float sumChargedHadronPt
sum-pt of charged Hadron // old float chargedHadronIso ;
Definition: GsfElectron.h:559

Member Data Documentation

const float GsfEleEffAreaPFIsoCut::_barrelCutOff
private

Definition at line 26 of file GsfEleEffAreaPFIsoCut.cc.

Referenced by operator()().

EffectiveAreas GsfEleEffAreaPFIsoCut::_effectiveAreas
private

Definition at line 29 of file GsfEleEffAreaPFIsoCut.cc.

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

const float GsfEleEffAreaPFIsoCut::_isoCutEBHighPt
private

Definition at line 23 of file GsfEleEffAreaPFIsoCut.cc.

Referenced by operator()().

const float GsfEleEffAreaPFIsoCut::_isoCutEBLowPt
private

Definition at line 23 of file GsfEleEffAreaPFIsoCut.cc.

Referenced by operator()().

const float GsfEleEffAreaPFIsoCut::_isoCutEEHighPt
private

Definition at line 23 of file GsfEleEffAreaPFIsoCut.cc.

Referenced by operator()().

const float GsfEleEffAreaPFIsoCut::_isoCutEELowPt
private

Definition at line 23 of file GsfEleEffAreaPFIsoCut.cc.

Referenced by operator()().

bool GsfEleEffAreaPFIsoCut::_isRelativeIso
private

Definition at line 27 of file GsfEleEffAreaPFIsoCut.cc.

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

const float GsfEleEffAreaPFIsoCut::_ptCutOff
private

Definition at line 25 of file GsfEleEffAreaPFIsoCut.cc.

Referenced by operator()().

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

Definition at line 31 of file GsfEleEffAreaPFIsoCut.cc.

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

char GsfEleEffAreaPFIsoCut::rhoString_ = "rho"
staticprivate

Definition at line 33 of file GsfEleEffAreaPFIsoCut.cc.

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