Main Page
Namespaces
Classes
Package Documentation
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
7
#include "
RecoEgamma/ElectronIdentification/interface/EBEECutValues.h
"
8
9
class
GsfEleTrkPtIsoRhoCut
:
public
CutApplicatorWithEventContentBase
{
10
public
:
11
GsfEleTrkPtIsoRhoCut
(
const
edm::ParameterSet
&
c
);
12
13
result_type
operator()
(
const
reco::GsfElectronPtr
&)
const
final
;
14
15
void
setConsumes
(
edm::ConsumesCollector
&)
final
;
16
void
getEventContent
(
const
edm::EventBase
&)
final
;
17
18
double
value
(
const
reco::CandidatePtr
&
cand
)
const
final
;
19
20
CandidateType
candidateType
() const final {
21
return
ELECTRON
;
22
}
23
24
private
:
25
26
EBEECutValues
slopeTerm_
;
27
EBEECutValues
slopeStart_
;
28
EBEECutValues
constTerm_
;
29
EBEECutValues
rhoEtStart_
;
30
EBEECutValues
rhoEA_
;
31
32
33
edm::Handle<double>
rhoHandle_
;
34
35
};
36
37
DEFINE_EDM_PLUGIN
(
CutApplicatorFactory
,
38
GsfEleTrkPtIsoRhoCut
,
39
"GsfEleTrkPtIsoRhoCut"
);
40
41
GsfEleTrkPtIsoRhoCut::GsfEleTrkPtIsoRhoCut
(
const
edm::ParameterSet
& params) :
42
CutApplicatorWithEventContentBase
(params),
43
slopeTerm_
(params,
"slopeTerm"
),
44
slopeStart_
(params,
"slopeStart"
),
45
constTerm_
(params,
"constTerm"
),
46
rhoEtStart_
(params,
"rhoEtStart"
),
47
rhoEA_
(params,
"rhoEA"
)
48
{
49
edm::InputTag
rhoTag
= params.
getParameter
<
edm::InputTag
>(
"rho"
);
50
contentTags_
.emplace(
"rho"
,rhoTag);
51
}
52
53
void
GsfEleTrkPtIsoRhoCut::setConsumes
(
edm::ConsumesCollector
& cc) {
54
auto
rho
= cc.
consumes
<
double
>(
contentTags_
[
"rho"
]);
55
contentTokens_
.emplace(
"rho"
,
rho
);
56
}
57
58
void
GsfEleTrkPtIsoRhoCut::getEventContent
(
const
edm::EventBase
&
ev
) {
59
ev.
getByLabel
(
contentTags_
[
"rho"
],
rhoHandle_
);
60
}
61
62
CutApplicatorBase::result_type
63
GsfEleTrkPtIsoRhoCut::
64
operator()
(
const
reco::GsfElectronPtr
&
cand
)
const
{
65
const
double
rho
= (*rhoHandle_);
66
const
float
isolTrkPt = cand->
dr03TkSumPt
();
67
68
const
float
et
= cand->
et
();
69
const
float
cutValue = et >
slopeStart_
(cand) ?
slopeTerm_
(cand)*(et-
slopeStart_
(cand)) +
constTerm_
(cand) :
constTerm_
(cand);
70
71
const
float
rhoCutValue = et >=
rhoEtStart_
(cand) ?
rhoEA_
(cand)*rho : 0.;
72
73
return
isolTrkPt < cutValue+rhoCutValue;
74
}
75
76
double
GsfEleTrkPtIsoRhoCut::
77
value
(
const
reco::CandidatePtr
&
cand
)
const
{
78
reco::GsfElectronPtr
ele(cand);
79
return
ele->
dr03TkSumPt
();
80
}
edm::ConsumesCollector::consumes
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
Definition:
ConsumesCollector.h:52
DDAxes::rho
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
GsfEleTrkPtIsoRhoCut::value
double value(const reco::CandidatePtr &cand) const final
Definition:
GsfEleTrkPtIsoRhoCut.cc:77
CutApplicatorWithEventContentBase
Definition:
CutApplicatorWithEventContentBase.h:19
CutApplicatorBase::ELECTRON
Definition:
CutApplicatorBase.h:48
CutApplicatorWithEventContentBase::contentTokens_
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
Definition:
CutApplicatorWithEventContentBase.h:44
GsfEleTrkPtIsoRhoCut::candidateType
CandidateType candidateType() const final
Definition:
GsfEleTrkPtIsoRhoCut.cc:20
EnergyCorrector.c
c
Definition:
EnergyCorrector.py:43
edm::Handle< double >
GsfEleTrkPtIsoRhoCut::constTerm_
EBEECutValues constTerm_
Definition:
GsfEleTrkPtIsoRhoCut.cc:28
Conversion.h
ev
bool ev
Definition:
Hydjet2Hadronizer.cc:95
EBEECutValuesT< double >
GsfEleTrkPtIsoRhoCut::rhoEA_
EBEECutValues rhoEA_
Definition:
GsfEleTrkPtIsoRhoCut.cc:30
CutApplicatorWithEventContentBase::contentTags_
std::unordered_map< std::string, edm::InputTag > contentTags_
Definition:
CutApplicatorWithEventContentBase.h:39
GsfEleTrkPtIsoRhoCut::getEventContent
void getEventContent(const edm::EventBase &) final
Definition:
GsfEleTrkPtIsoRhoCut.cc:58
edmplugin::PluginFactory
Definition:
PluginFactory.h:32
reco::LeafCandidate::et
double et() const final
transverse energy
Definition:
LeafCandidate.h:112
ConversionTools.h
GsfEleTrkPtIsoRhoCut::rhoHandle_
edm::Handle< double > rhoHandle_
Definition:
GsfEleTrkPtIsoRhoCut.cc:33
GsfEleTrkPtIsoRhoCut::setConsumes
void setConsumes(edm::ConsumesCollector &) final
Definition:
GsfEleTrkPtIsoRhoCut.cc:53
reco::GsfElectron::dr03TkSumPt
float dr03TkSumPt() const
Definition:
GsfElectron.h:543
CutApplicatorWithEventContentBase.h
edm::Ptr< reco::GsfElectron >
GsfEleTrkPtIsoRhoCut::slopeTerm_
EBEECutValues slopeTerm_
Definition:
GsfEleTrkPtIsoRhoCut.cc:26
GsfEleTrkPtIsoRhoCut::GsfEleTrkPtIsoRhoCut
GsfEleTrkPtIsoRhoCut(const edm::ParameterSet &c)
Definition:
GsfEleTrkPtIsoRhoCut.cc:41
edm::EventBase
Definition:
EventBase.h:46
GsfElectron.h
EBEECutValues.h
GsfEleTrkPtIsoRhoCut::operator()
result_type operator()(const reco::GsfElectronPtr &) const final
Definition:
GsfEleTrkPtIsoRhoCut.cc:64
stringResolutionProvider_cfi.et
et
define resolution functions of each parameter
Definition:
stringResolutionProvider_cfi.py:13
GsfEleTrkPtIsoRhoCut::slopeStart_
EBEECutValues slopeStart_
Definition:
GsfEleTrkPtIsoRhoCut.cc:27
ConversionFwd.h
edm::EventBase::getByLabel
bool getByLabel(InputTag const &, Handle< T > &) const
Definition:
EventBase.h:94
GsfEleTrkPtIsoRhoCut
Definition:
GsfEleTrkPtIsoRhoCut.cc:9
edm::InputTag
Definition:
InputTag.h:15
CutApplicatorBase::CandidateType
CandidateType
Definition:
CutApplicatorBase.h:47
GsfEleTrkPtIsoRhoCut::rhoEtStart_
EBEECutValues rhoEtStart_
Definition:
GsfEleTrkPtIsoRhoCut.cc:29
edm::ParameterSet
Definition:
ParameterSet.h:36
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition:
PluginFactory.h:102
cand
Definition:
decayParser.h:34
edm::ConsumesCollector
Definition:
ConsumesCollector.h:39
HiL1Corrector_cff.rhoTag
rhoTag
Definition:
HiL1Corrector_cff.py:7
Generated for CMSSW Reference Manual by
1.8.11