RecoEgamma
ElectronIdentification
plugins
cuts
GsfEleTrkPtIsoRhoCut.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/EgammaCandidates/interface/ConversionFwd.h
"
4
#include "
DataFormats/EgammaCandidates/interface/Conversion.h
"
5
#include "
RecoEgamma/EgammaTools/interface/ConversionTools.h
"
6
#include "
RecoEgamma/EgammaTools/interface/EBEECutValues.h
"
7
8
class
GsfEleTrkPtIsoRhoCut
:
public
CutApplicatorWithEventContentBase
{
9
public
:
10
GsfEleTrkPtIsoRhoCut
(
const
edm::ParameterSet
&
c
);
11
12
result_type
operator()
(
const
reco::GsfElectronPtr
&)
const
final
;
13
14
void
setConsumes
(
edm::ConsumesCollector
&)
final
;
15
void
getEventContent
(
const
edm::EventBase
&)
final
;
16
17
double
value
(
const
reco::CandidatePtr
&
cand
)
const
final
;
18
19
CandidateType
candidateType
()
const
final {
return
ELECTRON
; }
20
21
private
:
22
EBEECutValues
slopeTerm_
;
23
EBEECutValues
slopeStart_
;
24
EBEECutValues
constTerm_
;
25
EBEECutValues
rhoEtStart_
;
26
EBEECutValues
rhoEA_
;
27
28
edm::Handle<double>
rhoHandle_
;
29
};
30
31
DEFINE_EDM_PLUGIN
(
CutApplicatorFactory
,
GsfEleTrkPtIsoRhoCut
,
"GsfEleTrkPtIsoRhoCut"
);
32
33
GsfEleTrkPtIsoRhoCut::GsfEleTrkPtIsoRhoCut
(
const
edm::ParameterSet
&
params
)
34
:
CutApplicatorWithEventContentBase
(
params
),
35
slopeTerm_(
params
,
"slopeTerm"
),
36
slopeStart_(
params
,
"slopeStart"
),
37
constTerm_(
params
,
"constTerm"
),
38
rhoEtStart_(
params
,
"rhoEtStart"
),
39
rhoEA_(
params
,
"rhoEA"
) {
40
edm::InputTag
rhoTag
=
params
.getParameter<
edm::InputTag
>(
"rho"
);
41
contentTags_
.emplace(
"rho"
,
rhoTag
);
42
}
43
44
void
GsfEleTrkPtIsoRhoCut::setConsumes
(
edm::ConsumesCollector
&
cc
) {
45
auto
rho
=
cc
.consumes<
double
>(
contentTags_
[
"rho"
]);
46
contentTokens_
.emplace(
"rho"
,
rho
);
47
}
48
49
void
GsfEleTrkPtIsoRhoCut::getEventContent
(
const
edm::EventBase
&
ev
) {
ev
.getByLabel(
contentTags_
[
"rho"
],
rhoHandle_
); }
50
51
CutApplicatorBase::result_type
GsfEleTrkPtIsoRhoCut::operator()
(
const
reco::GsfElectronPtr
&
cand
)
const
{
52
const
double
rho
= (*rhoHandle_);
53
const
float
isolTrkPt =
cand
->dr03TkSumPt();
54
55
const
float
et
=
cand
->et();
56
const
float
cutValue =
57
et
>
slopeStart_
(
cand
) ?
slopeTerm_
(
cand
) * (
et
-
slopeStart_
(
cand
)) +
constTerm_
(
cand
) :
constTerm_
(
cand
);
58
59
const
float
rhoCutValue =
et
>=
rhoEtStart_
(
cand
) ?
rhoEA_
(
cand
) *
rho
: 0.;
60
61
return
isolTrkPt < cutValue + rhoCutValue;
62
}
63
64
double
GsfEleTrkPtIsoRhoCut::value
(
const
reco::CandidatePtr
&
cand
)
const
{
65
reco::GsfElectronPtr
ele(
cand
);
66
return
ele->
dr03TkSumPt
();
67
}
GsfEleTrkPtIsoRhoCut::candidateType
CandidateType candidateType() const final
Definition:
GsfEleTrkPtIsoRhoCut.cc:19
reco::GsfElectron::dr03TkSumPt
float dr03TkSumPt() const
Definition:
GsfElectron.h:529
ConversionTools.h
CalibrationSummaryClient_cfi.params
params
Definition:
CalibrationSummaryClient_cfi.py:14
EBEECutValues.h
HLT_2018_cff.rhoTag
rhoTag
Definition:
HLT_2018_cff.py:13606
ConversionFwd.h
watchdog.const
const
Definition:
watchdog.py:83
edm::Handle< double >
CutApplicatorWithEventContentBase::contentTags_
std::unordered_map< std::string, edm::InputTag > contentTags_
Definition:
CutApplicatorWithEventContentBase.h:35
candidate_functions::CandidateCut::result_type
bool result_type
Definition:
CandidateCut.h:11
GsfEleTrkPtIsoRhoCut::value
double value(const reco::CandidatePtr &cand) const final
Definition:
GsfEleTrkPtIsoRhoCut.cc:64
CutApplicatorWithEventContentBase.h
GsfElectron.h
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition:
PluginFactory.h:124
GsfEleTrkPtIsoRhoCut::getEventContent
void getEventContent(const edm::EventBase &) final
Definition:
GsfEleTrkPtIsoRhoCut.cc:49
GsfEleTrkPtIsoRhoCut::rhoEA_
EBEECutValues rhoEA_
Definition:
GsfEleTrkPtIsoRhoCut.cc:26
DDAxes::rho
GsfEleTrkPtIsoRhoCut::rhoHandle_
edm::Handle< double > rhoHandle_
Definition:
GsfEleTrkPtIsoRhoCut.cc:28
GsfEleTrkPtIsoRhoCut::GsfEleTrkPtIsoRhoCut
GsfEleTrkPtIsoRhoCut(const edm::ParameterSet &c)
Definition:
GsfEleTrkPtIsoRhoCut.cc:33
edm::ParameterSet
Definition:
ParameterSet.h:36
CutApplicatorBase::CandidateType
CandidateType
Definition:
CutApplicatorBase.h:47
GsfEleTrkPtIsoRhoCut::constTerm_
EBEECutValues constTerm_
Definition:
GsfEleTrkPtIsoRhoCut.cc:24
edmplugin::PluginFactory
Definition:
PluginFactory.h:34
GsfEleTrkPtIsoRhoCut::operator()
result_type operator()(const reco::GsfElectronPtr &) const final
Definition:
GsfEleTrkPtIsoRhoCut.cc:51
cand
Definition:
decayParser.h:34
EgHLTOffHistBins_cfi.et
et
Definition:
EgHLTOffHistBins_cfi.py:8
HltBtagPostValidation_cff.c
c
Definition:
HltBtagPostValidation_cff.py:31
cc
edm::Ptr< reco::GsfElectron >
GsfEleTrkPtIsoRhoCut::slopeStart_
EBEECutValues slopeStart_
Definition:
GsfEleTrkPtIsoRhoCut.cc:23
GsfEleTrkPtIsoRhoCut::slopeTerm_
EBEECutValues slopeTerm_
Definition:
GsfEleTrkPtIsoRhoCut.cc:22
GsfEleTrkPtIsoRhoCut::rhoEtStart_
EBEECutValues rhoEtStart_
Definition:
GsfEleTrkPtIsoRhoCut.cc:25
GsfEleTrkPtIsoRhoCut::setConsumes
void setConsumes(edm::ConsumesCollector &) final
Definition:
GsfEleTrkPtIsoRhoCut.cc:44
ev
bool ev
Definition:
Hydjet2Hadronizer.cc:95
CutApplicatorWithEventContentBase::contentTokens_
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
Definition:
CutApplicatorWithEventContentBase.h:40
edm::EventBase
Definition:
EventBase.h:46
EBEECutValuesT< double >
CutApplicatorBase::ELECTRON
Definition:
CutApplicatorBase.h:47
edm::InputTag
Definition:
InputTag.h:15
edm::ConsumesCollector
Definition:
ConsumesCollector.h:39
GsfEleTrkPtIsoRhoCut
Definition:
GsfEleTrkPtIsoRhoCut.cc:8
CutApplicatorWithEventContentBase
Definition:
CutApplicatorWithEventContentBase.h:19
Conversion.h
Generated for CMSSW Reference Manual by
1.8.16