src
RecoEgamma
ElectronIdentification
plugins
cuts
GsfEleHadronicOverEMEnergyScaledCut.cc
Go to the documentation of this file.
1
#include "
PhysicsTools/SelectorUtils/interface/CutApplicatorWithEventContentBase.h
"
2
#include "
DataFormats/EgammaCandidates/interface/GsfElectron.h
"
3
4
class
GsfEleHadronicOverEMEnergyScaledCut
:
public
CutApplicatorWithEventContentBase
{
5
public
:
6
GsfEleHadronicOverEMEnergyScaledCut
(
const
edm::ParameterSet
&
c
)
7
:
CutApplicatorWithEventContentBase
(
c
),
8
barrelC0_
(
c
.getParameter<double>(
"barrelC0"
)),
9
barrelCE_
(
c
.getParameter<double>(
"barrelCE"
)),
10
barrelCr_
(
c
.getParameter<double>(
"barrelCr"
)),
11
endcapC0_
(
c
.getParameter<double>(
"endcapC0"
)),
12
endcapCE_
(
c
.getParameter<double>(
"endcapCE"
)),
13
endcapCr_
(
c
.getParameter<double>(
"endcapCr"
)),
14
barrelCutOff_
(
c
.getParameter<double>(
"barrelCutOff"
)) {
15
edm::InputTag
rhoTag
=
c
.getParameter<
edm::InputTag
>(
"rho"
);
16
contentTags_
.emplace(
"rho"
,
rhoTag
);
17
}
18
19
result_type
operator()
(
const
reco::GsfElectronPtr
&)
const
final
;
20
21
void
setConsumes
(
edm::ConsumesCollector
&)
final
;
22
void
getEventContent
(
const
edm::EventBase
&)
final
;
23
24
double
value
(
const
reco::CandidatePtr
&
cand
)
const
final
;
25
26
CandidateType
candidateType
()
const
final {
return
ELECTRON
; }
27
28
private
:
29
const
float
barrelC0_
,
barrelCE_
,
barrelCr_
,
endcapC0_
,
endcapCE_
,
endcapCr_
,
barrelCutOff_
;
30
edm::Handle<double>
rhoHandle_
;
31
};
32
33
DEFINE_EDM_PLUGIN
(
CutApplicatorFactory
,
GsfEleHadronicOverEMEnergyScaledCut
,
"GsfEleHadronicOverEMEnergyScaledCut"
);
34
35
void
GsfEleHadronicOverEMEnergyScaledCut::setConsumes
(
edm::ConsumesCollector
&
cc
) {
36
auto
rho
=
cc
.consumes<
double
>(
contentTags_
[
"rho"
]);
37
contentTokens_
.emplace(
"rho"
,
rho
);
38
}
39
40
void
GsfEleHadronicOverEMEnergyScaledCut::getEventContent
(
const
edm::EventBase
&
ev
) {
41
ev
.getByLabel(
contentTags_
[
"rho"
],
rhoHandle_
);
42
}
43
44
CutApplicatorBase::result_type
GsfEleHadronicOverEMEnergyScaledCut::operator()
(
const
reco::GsfElectronPtr
&
cand
)
const
{
45
const
double
rho
=
rhoHandle_
.
isValid
() ? (
float
)(*
rhoHandle_
) : 0;
46
const
float
energy
=
cand
->superCluster()->energy();
47
const
float
c0
= (
std::abs
(
cand
->superCluster()->position().eta()) <
barrelCutOff_
?
barrelC0_
:
endcapC0_
);
48
const
float
cE = (
std::abs
(
cand
->superCluster()->position().eta()) <
barrelCutOff_
?
barrelCE_
:
endcapCE_
);
49
const
float
cR = (
std::abs
(
cand
->superCluster()->position().eta()) <
barrelCutOff_
?
barrelCr_
:
endcapCr_
);
50
return
cand
->hadronicOverEm() <
c0
+ cE /
energy
+ cR *
rho
/
energy
;
51
}
52
53
double
GsfEleHadronicOverEMEnergyScaledCut::value
(
const
reco::CandidatePtr
&
cand
)
const
{
54
reco::GsfElectronPtr
ele(
cand
);
55
return
ele->
hadronicOverEm
();
56
}
DDAxes::rho
makeMEIFBenchmarkPlots.ev
ev
Definition:
makeMEIFBenchmarkPlots.py:55
CutApplicatorWithEventContentBase
Definition:
CutApplicatorWithEventContentBase.h:19
HLT_2024v12_cff.rhoTag
rhoTag
Definition:
HLT_2024v12_cff.py:17874
GsfEleHadronicOverEMEnergyScaledCut::operator()
result_type operator()(const reco::GsfElectronPtr &) const final
Definition:
GsfEleHadronicOverEMEnergyScaledCut.cc:44
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:374
GsfEleHadronicOverEMEnergyScaledCut::value
double value(const reco::CandidatePtr &cand) const final
Definition:
GsfEleHadronicOverEMEnergyScaledCut.cc:53
gpuPixelDoublets::cc
uint32_t cc[maxCellsPerHit]
Definition:
gpuFishbone.h:49
edm::Handle< double >
GsfEleHadronicOverEMEnergyScaledCut::barrelCE_
const float barrelCE_
Definition:
GsfEleHadronicOverEMEnergyScaledCut.cc:29
GsfEleHadronicOverEMEnergyScaledCut::barrelCr_
const float barrelCr_
Definition:
GsfEleHadronicOverEMEnergyScaledCut.cc:29
GsfEleHadronicOverEMEnergyScaledCut::barrelC0_
const float barrelC0_
Definition:
GsfEleHadronicOverEMEnergyScaledCut.cc:29
hcalRecHitTable_cff.energy
energy
Definition:
hcalRecHitTable_cff.py:13
HltBtagPostValidation_cff.c
c
Definition:
HltBtagPostValidation_cff.py:35
GsfEleHadronicOverEMEnergyScaledCut::barrelCutOff_
const float barrelCutOff_
Definition:
GsfEleHadronicOverEMEnergyScaledCut.cc:29
nano_mu_digi_cff.float
float
Definition:
nano_mu_digi_cff.py:13
GsfEleHadronicOverEMEnergyScaledCut::GsfEleHadronicOverEMEnergyScaledCut
GsfEleHadronicOverEMEnergyScaledCut(const edm::ParameterSet &c)
Definition:
GsfEleHadronicOverEMEnergyScaledCut.cc:6
edmplugin::PluginFactory
Definition:
PluginFactory.h:35
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
CutApplicatorWithEventContentBase.h
edm::Ptr< reco::GsfElectron >
GsfEleHadronicOverEMEnergyScaledCut
Definition:
GsfEleHadronicOverEMEnergyScaledCut.cc:4
edm::EventBase
Definition:
EventBase.h:47
GsfEleHadronicOverEMEnergyScaledCut::endcapC0_
const float endcapC0_
Definition:
GsfEleHadronicOverEMEnergyScaledCut.cc:29
GsfElectron.h
GsfEleHadronicOverEMEnergyScaledCut::candidateType
CandidateType candidateType() const final
Definition:
GsfEleHadronicOverEMEnergyScaledCut.cc:26
CutApplicatorWithEventContentBase::contentTokens_
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
Definition:
CutApplicatorWithEventContentBase.h:40
GsfEleHadronicOverEMEnergyScaledCut::setConsumes
void setConsumes(edm::ConsumesCollector &) final
Definition:
GsfEleHadronicOverEMEnergyScaledCut.cc:35
edm::HandleBase::isValid
bool isValid() const
Definition:
HandleBase.h:70
edm::InputTag
Definition:
InputTag.h:15
GsfEleHadronicOverEMEnergyScaledCut::endcapCE_
const float endcapCE_
Definition:
GsfEleHadronicOverEMEnergyScaledCut.cc:29
CutApplicatorBase::CandidateType
CandidateType
Definition:
CutApplicatorBase.h:47
GsfEleHadronicOverEMEnergyScaledCut::rhoHandle_
edm::Handle< double > rhoHandle_
Definition:
GsfEleHadronicOverEMEnergyScaledCut.cc:30
fftjetpileupestimator_calo_uncalib_cfi.c0
c0
Definition:
fftjetpileupestimator_calo_uncalib_cfi.py:8
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
reco::GsfElectron::hadronicOverEm
float hadronicOverEm() const
Definition:
GsfElectron.h:500
GsfEleHadronicOverEMEnergyScaledCut::getEventContent
void getEventContent(const edm::EventBase &) final
Definition:
GsfEleHadronicOverEMEnergyScaledCut.cc:40
GsfEleHadronicOverEMEnergyScaledCut::endcapCr_
const float endcapCr_
Definition:
GsfEleHadronicOverEMEnergyScaledCut.cc:29
edm::ConsumesCollector
Definition:
ConsumesCollector.h:45
Generated for CMSSW Reference Manual by
1.8.14