RecoEgamma
ElectronIdentification
plugins
cuts
GsfEleHadronicOverEMLinearCut.cc
Go to the documentation of this file.
1
#include "
PhysicsTools/SelectorUtils/interface/CutApplicatorBase.h
"
2
#include "
DataFormats/EgammaCandidates/interface/GsfElectron.h
"
3
#include "
RecoEgamma/EgammaTools/interface/EBEECutValues.h
"
4
5
class
GsfEleHadronicOverEMLinearCut
:
public
CutApplicatorBase
{
6
public
:
7
GsfEleHadronicOverEMLinearCut
(
const
edm::ParameterSet
&
params
)
8
:
CutApplicatorBase
(
params
),
9
slopeTerm_
(
params
,
"slopeTerm"
),
10
slopeStart_
(
params
,
"slopeStart"
),
11
constTerm_
(
params
,
"constTerm"
) {}
12
13
result_type
operator()
(
const
reco::GsfElectronPtr
&)
const
final
;
14
15
double
value
(
const
reco::CandidatePtr
&
cand
)
const
final
;
16
17
CandidateType
candidateType
()
const
final {
return
ELECTRON
; }
18
19
private
:
20
EBEECutValues
slopeTerm_
;
21
EBEECutValues
slopeStart_
;
22
EBEECutValues
constTerm_
;
23
};
24
25
DEFINE_EDM_PLUGIN
(
CutApplicatorFactory
,
GsfEleHadronicOverEMLinearCut
,
"GsfEleHadronicOverEMLinearCut"
);
26
27
CutApplicatorBase::result_type
GsfEleHadronicOverEMLinearCut::operator()
(
const
reco::GsfElectronPtr
&
cand
)
const
{
28
const
float
energy
=
cand
->superCluster()->energy();
29
const
float
cutValue =
energy
>
slopeStart_
(
cand
) ?
slopeTerm_
(
cand
) * (
energy
-
slopeStart_
(
cand
)) +
constTerm_
(
cand
)
30
:
constTerm_
(
cand
);
31
32
return
cand
->hadronicOverEm() *
energy
< cutValue;
33
}
34
35
double
GsfEleHadronicOverEMLinearCut::value
(
const
reco::CandidatePtr
&
cand
)
const
{
36
reco::GsfElectronPtr
ele(
cand
);
37
const
float
energy
= ele->
superCluster
()->energy();
38
return
ele->
hadronicOverEm
() *
energy
;
39
}
CalibrationSummaryClient_cfi.params
params
Definition:
CalibrationSummaryClient_cfi.py:14
EBEECutValues.h
GsfEleHadronicOverEMLinearCut::slopeTerm_
EBEECutValues slopeTerm_
Definition:
GsfEleHadronicOverEMLinearCut.cc:20
reco::GsfElectron::hadronicOverEm
float hadronicOverEm() const
Definition:
GsfElectron.h:508
watchdog.const
const
Definition:
watchdog.py:83
candidate_functions::CandidateCut::result_type
bool result_type
Definition:
CandidateCut.h:11
GsfEleHadronicOverEMLinearCut::value
double value(const reco::CandidatePtr &cand) const final
Definition:
GsfEleHadronicOverEMLinearCut.cc:35
HCALHighEnergyHPDFilter_cfi.energy
energy
Definition:
HCALHighEnergyHPDFilter_cfi.py:5
GsfEleHadronicOverEMLinearCut::candidateType
CandidateType candidateType() const final
Definition:
GsfEleHadronicOverEMLinearCut.cc:17
GsfElectron.h
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition:
PluginFactory.h:124
edm::ParameterSet
Definition:
ParameterSet.h:47
GsfEleHadronicOverEMLinearCut::GsfEleHadronicOverEMLinearCut
GsfEleHadronicOverEMLinearCut(const edm::ParameterSet ¶ms)
Definition:
GsfEleHadronicOverEMLinearCut.cc:7
GsfEleHadronicOverEMLinearCut::operator()
result_type operator()(const reco::GsfElectronPtr &) const final
Definition:
GsfEleHadronicOverEMLinearCut.cc:27
CutApplicatorBase::CandidateType
CandidateType
Definition:
CutApplicatorBase.h:47
edmplugin::PluginFactory
Definition:
PluginFactory.h:34
cand
Definition:
decayParser.h:32
GsfEleHadronicOverEMLinearCut::constTerm_
EBEECutValues constTerm_
Definition:
GsfEleHadronicOverEMLinearCut.cc:22
GsfEleHadronicOverEMLinearCut::slopeStart_
EBEECutValues slopeStart_
Definition:
GsfEleHadronicOverEMLinearCut.cc:21
edm::Ptr< reco::GsfElectron >
CutApplicatorBase
Definition:
CutApplicatorBase.h:45
GsfEleHadronicOverEMLinearCut
Definition:
GsfEleHadronicOverEMLinearCut.cc:5
EBEECutValuesT< double >
reco::GsfElectron::superCluster
SuperClusterRef superCluster() const override
reference to a SuperCluster
Definition:
GsfElectron.h:163
CutApplicatorBase.h
CutApplicatorBase::ELECTRON
Definition:
CutApplicatorBase.h:47
Generated for CMSSW Reference Manual by
1.8.16