src
RecoEgamma
ElectronIdentification
plugins
cuts
GsfEleRelPFIsoScaledCut.cc
Go to the documentation of this file.
1
#include "
PhysicsTools/SelectorUtils/interface/CutApplicatorWithEventContentBase.h
"
2
#include "
DataFormats/EgammaCandidates/interface/GsfElectron.h
"
3
#include "
CommonTools/Egamma/interface/EffectiveAreas.h
"
4
5
class
GsfEleRelPFIsoScaledCut
:
public
CutApplicatorWithEventContentBase
{
6
public
:
7
GsfEleRelPFIsoScaledCut
(
const
edm::ParameterSet
&
c
);
8
9
result_type
operator()
(
const
reco::GsfElectronPtr
&)
const
final
;
10
11
void
setConsumes
(
edm::ConsumesCollector
&)
final
;
12
void
getEventContent
(
const
edm::EventBase
&)
final
;
13
14
double
value
(
const
reco::CandidatePtr
&
cand
)
const
final
;
15
16
CandidateType
candidateType
()
const
final {
return
ELECTRON
; }
17
18
private
:
19
const
float
barrelC0_
,
endcapC0_
,
barrelCpt_
,
endcapCpt_
,
barrelCutOff_
;
20
EffectiveAreas
effectiveAreas_
;
21
edm::Handle<double>
rhoHandle_
;
22
};
23
24
DEFINE_EDM_PLUGIN
(
CutApplicatorFactory
,
GsfEleRelPFIsoScaledCut
,
"GsfEleRelPFIsoScaledCut"
);
25
26
GsfEleRelPFIsoScaledCut::GsfEleRelPFIsoScaledCut
(
const
edm::ParameterSet
&
c
)
27
:
CutApplicatorWithEventContentBase
(
c
),
28
barrelC0_(
c
.getParameter<double>(
"barrelC0"
)),
29
endcapC0_(
c
.getParameter<double>(
"endcapC0"
)),
30
barrelCpt_(
c
.getParameter<double>(
"barrelCpt"
)),
31
endcapCpt_(
c
.getParameter<double>(
"endcapCpt"
)),
32
barrelCutOff_(
c
.getParameter<double>(
"barrelCutOff"
)),
33
effectiveAreas_((
c
.getParameter<
edm
::FileInPath>(
"effAreasConfigFile"
)).
fullPath
()) {
34
edm::InputTag
rhoTag
=
c
.getParameter<
edm::InputTag
>(
"rho"
);
35
contentTags_
.emplace(
"rho"
,
rhoTag
);
36
}
37
38
void
GsfEleRelPFIsoScaledCut::setConsumes
(
edm::ConsumesCollector
&
cc
) {
39
auto
rho
=
cc
.consumes<
double
>(
contentTags_
[
"rho"
]);
40
contentTokens_
.emplace(
"rho"
,
rho
);
41
}
42
43
void
GsfEleRelPFIsoScaledCut::getEventContent
(
const
edm::EventBase
&
ev
) {
44
ev
.getByLabel(
contentTags_
[
"rho"
],
rhoHandle_
);
45
}
46
47
CutApplicatorBase::result_type
GsfEleRelPFIsoScaledCut::operator()
(
const
reco::GsfElectronPtr
&
cand
)
const
{
48
// Establish the cut value
49
double
absEta =
std::abs
(
cand
->superCluster()->eta());
50
51
const
float
C0 = (absEta <
barrelCutOff_
?
barrelC0_
:
endcapC0_
);
52
const
float
Cpt = (absEta <
barrelCutOff_
?
barrelCpt_
:
endcapCpt_
);
53
const
float
isoCut = C0 + Cpt /
cand
->pt();
54
55
return
value
(
cand
) < isoCut;
56
}
57
58
double
GsfEleRelPFIsoScaledCut::value
(
const
reco::CandidatePtr
&
cand
)
const
{
59
// Establish the cut value
60
reco::GsfElectronPtr
ele(
cand
);
61
double
absEta =
std::abs
(ele->
superCluster
()->eta());
62
63
// Compute the combined isolation with effective area correction
64
auto
pfIso = ele->
pfIsolationVariables
();
65
const
float
chad = pfIso.
sumChargedHadronPt
;
66
const
float
nhad = pfIso.sumNeutralHadronEt;
67
const
float
pho = pfIso.sumPhotonEt;
68
const
float
eA =
effectiveAreas_
.
getEffectiveArea
(absEta);
69
const
float
rho
=
rhoHandle_
.
isValid
() ? (
float
)(*
rhoHandle_
) : 0;
// std::max likes float arguments
70
const
float
iso = chad +
std::max
(0.0
f
, nhad + pho -
rho
* eA);
71
return
iso /
cand
->pt();
72
}
DDAxes::rho
makeMEIFBenchmarkPlots.ev
ev
Definition:
makeMEIFBenchmarkPlots.py:55
CutApplicatorWithEventContentBase
Definition:
CutApplicatorWithEventContentBase.h:19
GsfEleRelPFIsoScaledCut::barrelCutOff_
const float barrelCutOff_
Definition:
GsfEleRelPFIsoScaledCut.cc:19
CutApplicatorWithEventContentBase::contentTags_
std::unordered_map< std::string, edm::InputTag > contentTags_
Definition:
CutApplicatorWithEventContentBase.h:35
CutApplicatorBase::ELECTRON
Definition:
CutApplicatorBase.h:47
runTheMatrix.const
const
Definition:
runTheMatrix.py:341
GsfEleRelPFIsoScaledCut::effectiveAreas_
EffectiveAreas effectiveAreas_
Definition:
GsfEleRelPFIsoScaledCut.cc:20
gpuPixelDoublets::cc
uint32_t cc[maxCellsPerHit]
Definition:
gpuFishbone.h:49
edm::Handle< double >
EffectiveAreas.h
HltBtagPostValidation_cff.c
c
Definition:
HltBtagPostValidation_cff.py:31
GsfEleRelPFIsoScaledCut::GsfEleRelPFIsoScaledCut
GsfEleRelPFIsoScaledCut(const edm::ParameterSet &c)
Definition:
GsfEleRelPFIsoScaledCut.cc:26
reco::GsfElectron::pfIsolationVariables
const PflowIsolationVariables & pfIsolationVariables() const
Definition:
GsfElectron.h:729
GsfEleRelPFIsoScaledCut::candidateType
CandidateType candidateType() const final
Definition:
GsfEleRelPFIsoScaledCut.cc:16
edmplugin::PluginFactory
Definition:
PluginFactory.h:35
GsfEleRelPFIsoScaledCut::barrelC0_
const float barrelC0_
Definition:
GsfEleRelPFIsoScaledCut.cc:19
candidate_functions::CandidateCut::result_type
bool result_type
Definition:
CandidateCut.h:11
funct::abs
Abs< T >::type abs(const T &t)
Definition:
Abs.h:22
f
double f[11][100]
Definition:
MuScleFitUtils.cc:78
HLT_2023v12_cff.rhoTag
rhoTag
Definition:
HLT_2023v12_cff.py:17938
GsfEleRelPFIsoScaledCut::rhoHandle_
edm::Handle< double > rhoHandle_
Definition:
GsfEleRelPFIsoScaledCut.cc:21
CutApplicatorWithEventContentBase.h
edm::Ptr< reco::GsfElectron >
SiStripPI::max
Definition:
SiStripPayloadInspectorHelper.h:178
GsfEleRelPFIsoScaledCut
Definition:
GsfEleRelPFIsoScaledCut.cc:5
edm::EventBase
Definition:
EventBase.h:47
EffectiveAreas::getEffectiveArea
const float getEffectiveArea(float eta) const
Definition:
EffectiveAreas.cc:76
GsfElectron.h
GsfEleRelPFIsoScaledCut::endcapCpt_
const float endcapCpt_
Definition:
GsfEleRelPFIsoScaledCut.cc:19
CutApplicatorWithEventContentBase::contentTokens_
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
Definition:
CutApplicatorWithEventContentBase.h:40
GsfEleRelPFIsoScaledCut::operator()
result_type operator()(const reco::GsfElectronPtr &) const final
Definition:
GsfEleRelPFIsoScaledCut.cc:47
GsfEleRelPFIsoScaledCut::setConsumes
void setConsumes(edm::ConsumesCollector &) final
Definition:
GsfEleRelPFIsoScaledCut.cc:38
GsfEleRelPFIsoScaledCut::endcapC0_
const float endcapC0_
Definition:
GsfEleRelPFIsoScaledCut.cc:19
edm::HandleBase::isValid
bool isValid() const
Definition:
HandleBase.h:70
edm
HLT enums.
Definition:
AlignableModifier.h:19
edm::InputTag
Definition:
InputTag.h:15
CutApplicatorBase::CandidateType
CandidateType
Definition:
CutApplicatorBase.h:47
GsfEleRelPFIsoScaledCut::getEventContent
void getEventContent(const edm::EventBase &) final
Definition:
GsfEleRelPFIsoScaledCut.cc:43
dqmMemoryStats.float
float
Definition:
dqmMemoryStats.py:127
edm::ParameterSet
Definition:
ParameterSet.h:47
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition:
PluginFactory.h:123
cand
Definition:
decayParser.h:32
GsfEleRelPFIsoScaledCut::barrelCpt_
const float barrelCpt_
Definition:
GsfEleRelPFIsoScaledCut.cc:19
GsfEleRelPFIsoScaledCut::value
double value(const reco::CandidatePtr &cand) const final
Definition:
GsfEleRelPFIsoScaledCut.cc:58
reco::GsfElectron::PflowIsolationVariables::sumChargedHadronPt
float sumChargedHadronPt
sum-pt of charged Hadron // old float chargedHadronIso ;
Definition:
GsfElectron.h:663
EffectiveAreas
Definition:
EffectiveAreas.h:7
reco::GsfElectron::superCluster
SuperClusterRef superCluster() const override
reference to a SuperCluster
Definition:
GsfElectron.h:155
edm::ConsumesCollector
Definition:
ConsumesCollector.h:45
contentValuesFiles.fullPath
fullPath
Definition:
contentValuesFiles.py:64
Generated for CMSSW Reference Manual by
1.8.14