src
RecoEgamma
ElectronIdentification
plugins
cuts
GsfEleValueMapIsoRhoCut.cc
Go to the documentation of this file.
1
#include "
PhysicsTools/SelectorUtils/interface/CutApplicatorWithEventContentBase.h
"
2
#include "
DataFormats/EgammaCandidates/interface/GsfElectron.h
"
3
#include "
DataFormats/Common/interface/ValueMap.h
"
4
#include "
RecoEgamma/EgammaTools/interface/EBEECutValues.h
"
5
6
class
GsfEleValueMapIsoRhoCut
:
public
CutApplicatorWithEventContentBase
{
7
public
:
8
GsfEleValueMapIsoRhoCut
(
const
edm::ParameterSet
&
c
);
9
10
result_type
operator()
(
const
reco::GsfElectronPtr
&)
const
final
;
11
12
void
setConsumes
(
edm::ConsumesCollector
&)
final
;
13
void
getEventContent
(
const
edm::EventBase
&)
final
;
14
15
double
value
(
const
reco::CandidatePtr
&
cand
)
const
final
;
16
17
CandidateType
candidateType
()
const
final {
return
ELECTRON
; }
18
19
private
:
20
float
getRhoCorr
(
const
reco::GsfElectronPtr
&
cand
)
const
;
21
22
EBEECutValues
slopeTerm_
;
23
EBEECutValues
slopeStart_
;
24
EBEECutValues
constTerm_
;
25
EBEECutValues
rhoEtStart_
;
26
EBEECutValues
rhoEA_
;
27
28
bool
useRho_
;
29
30
edm::Handle<double>
rhoHandle_
;
31
edm::Handle<edm::ValueMap<float>
>
valueHandle_
;
32
};
33
34
DEFINE_EDM_PLUGIN
(
CutApplicatorFactory
,
GsfEleValueMapIsoRhoCut
,
"GsfEleValueMapIsoRhoCut"
);
35
36
GsfEleValueMapIsoRhoCut::GsfEleValueMapIsoRhoCut
(
const
edm::ParameterSet
&
params
)
37
:
CutApplicatorWithEventContentBase
(
params
),
38
slopeTerm_(
params
,
"slopeTerm"
),
39
slopeStart_(
params
,
"slopeStart"
),
40
constTerm_(
params
,
"constTerm"
),
41
rhoEtStart_(
params
,
"rhoEtStart"
),
42
rhoEA_(
params
,
"rhoEA"
) {
43
auto
rho
=
params
.getParameter<
edm::InputTag
>(
"rho"
);
44
if
(!
rho
.label().empty()) {
45
useRho_
=
true
;
46
contentTags_
.emplace(
"rho"
,
rho
);
47
}
else
48
useRho_
=
false
;
49
50
contentTags_
.emplace(
"value"
,
params
.getParameter<
edm::InputTag
>(
"value"
));
51
}
52
53
void
GsfEleValueMapIsoRhoCut::setConsumes
(
edm::ConsumesCollector
&
cc
) {
54
if
(
useRho_
)
55
contentTokens_
.emplace(
"rho"
,
cc
.consumes<
double
>(
contentTags_
[
"rho"
]));
56
contentTokens_
.emplace(
"value"
,
cc
.consumes<
edm::ValueMap<float>
>(
contentTags_
[
"value"
]));
57
}
58
59
void
GsfEleValueMapIsoRhoCut::getEventContent
(
const
edm::EventBase
&
ev
) {
60
if
(
useRho_
)
61
ev
.getByLabel(
contentTags_
[
"rho"
],
rhoHandle_
);
62
ev
.getByLabel(
contentTags_
[
"value"
],
valueHandle_
);
63
}
64
65
CutApplicatorBase::result_type
GsfEleValueMapIsoRhoCut::operator()
(
const
reco::GsfElectronPtr
&
cand
)
const
{
66
const
float
val
= (*valueHandle_)[
cand
];
67
68
const
float
et
=
cand
->et();
69
const
float
cutValue =
70
et
>
slopeStart_
(
cand
) ?
slopeTerm_
(
cand
) * (
et
-
slopeStart_
(
cand
)) +
constTerm_
(
cand
) :
constTerm_
(
cand
);
71
const
float
rhoCutValue =
getRhoCorr
(
cand
);
72
73
return
val
< cutValue + rhoCutValue;
74
}
75
76
double
GsfEleValueMapIsoRhoCut::value
(
const
reco::CandidatePtr
&
cand
)
const
{
77
reco::GsfElectronPtr
ele(
cand
);
78
return
(*
valueHandle_
)[
cand
];
79
}
80
81
float
GsfEleValueMapIsoRhoCut::getRhoCorr
(
const
reco::GsfElectronPtr
&
cand
)
const
{
82
if
(!
useRho_
)
83
return
0.;
84
else
{
85
const
double
rho
= (*rhoHandle_);
86
return
cand
->et() >=
rhoEtStart_
(
cand
) ?
rhoEA_
(
cand
) *
rho
: 0.;
87
}
88
}
DDAxes::rho
makeMEIFBenchmarkPlots.ev
ev
Definition:
makeMEIFBenchmarkPlots.py:55
CutApplicatorWithEventContentBase
Definition:
CutApplicatorWithEventContentBase.h:19
GsfEleValueMapIsoRhoCut::setConsumes
void setConsumes(edm::ConsumesCollector &) final
Definition:
GsfEleValueMapIsoRhoCut.cc:53
DummyCfis.c
c
Definition:
DummyCfis.py:86
GsfEleValueMapIsoRhoCut::valueHandle_
edm::Handle< edm::ValueMap< float > > valueHandle_
Definition:
GsfEleValueMapIsoRhoCut.cc:31
CutApplicatorWithEventContentBase::contentTags_
std::unordered_map< std::string, edm::InputTag > contentTags_
Definition:
CutApplicatorWithEventContentBase.h:39
CutApplicatorBase::ELECTRON
Definition:
CutApplicatorBase.h:47
GsfEleValueMapIsoRhoCut::operator()
result_type operator()(const reco::GsfElectronPtr &) const final
Definition:
GsfEleValueMapIsoRhoCut.cc:65
gpuPixelDoublets::cc
uint32_t cc[maxCellsPerHit]
Definition:
gpuFishbone.h:49
GsfEleValueMapIsoRhoCut::value
double value(const reco::CandidatePtr &cand) const final
Definition:
GsfEleValueMapIsoRhoCut.cc:76
GsfEleValueMapIsoRhoCut::rhoEA_
EBEECutValues rhoEA_
Definition:
GsfEleValueMapIsoRhoCut.cc:26
edm::Handle< double >
EBEECutValuesT< double >
ValueMap.h
edmplugin::PluginFactory
Definition:
PluginFactory.h:35
GsfEleValueMapIsoRhoCut
Definition:
GsfEleValueMapIsoRhoCut.cc:6
GsfEleValueMapIsoRhoCut::useRho_
bool useRho_
Definition:
GsfEleValueMapIsoRhoCut.cc:28
candidate_functions::CandidateCut::result_type
bool result_type
Definition:
CandidateCut.h:11
GsfEleValueMapIsoRhoCut::getRhoCorr
float getRhoCorr(const reco::GsfElectronPtr &cand) const
Definition:
GsfEleValueMapIsoRhoCut.cc:81
GsfEleValueMapIsoRhoCut::slopeStart_
EBEECutValues slopeStart_
Definition:
GsfEleValueMapIsoRhoCut.cc:23
CutApplicatorWithEventContentBase.h
edm::Ptr< reco::GsfElectron >
edm::EventBase
Definition:
EventBase.h:47
submitPVValidationJobs.params
def params
Definition:
submitPVValidationJobs.py:483
GsfEleValueMapIsoRhoCut::rhoHandle_
edm::Handle< double > rhoHandle_
Definition:
GsfEleValueMapIsoRhoCut.cc:30
GsfElectron.h
GsfEleValueMapIsoRhoCut::rhoEtStart_
EBEECutValues rhoEtStart_
Definition:
GsfEleValueMapIsoRhoCut.cc:25
edm::ValueMap< float >
EBEECutValues.h
CutApplicatorWithEventContentBase::contentTokens_
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
Definition:
CutApplicatorWithEventContentBase.h:40
GsfEleValueMapIsoRhoCut::getEventContent
void getEventContent(const edm::EventBase &) final
Definition:
GsfEleValueMapIsoRhoCut.cc:59
Options.const
const
Definition:
Options.py:396
edm::InputTag
Definition:
InputTag.h:15
CutApplicatorBase::CandidateType
CandidateType
Definition:
CutApplicatorBase.h:47
GsfEleValueMapIsoRhoCut::candidateType
CandidateType candidateType() const final
Definition:
GsfEleValueMapIsoRhoCut.cc:17
edm::ParameterSet
Definition:
ParameterSet.h:48
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition:
PluginFactory.h:123
cand
Definition:
decayParser.h:32
heppy_batch.val
val
Definition:
heppy_batch.py:351
GsfEleValueMapIsoRhoCut::constTerm_
EBEECutValues constTerm_
Definition:
GsfEleValueMapIsoRhoCut.cc:24
GsfEleValueMapIsoRhoCut::GsfEleValueMapIsoRhoCut
GsfEleValueMapIsoRhoCut(const edm::ParameterSet &c)
Definition:
GsfEleValueMapIsoRhoCut.cc:36
l1tnanotables_cff.et
et
Definition:
l1tnanotables_cff.py:32
GsfEleValueMapIsoRhoCut::slopeTerm_
EBEECutValues slopeTerm_
Definition:
GsfEleValueMapIsoRhoCut.cc:22
edm::ConsumesCollector
Definition:
ConsumesCollector.h:45
Generated for CMSSW Reference Manual by
1.8.14