RecoEgamma
ElectronIdentification
plugins
cuts
GsfEleMinEcalEtCut.cc
Go to the documentation of this file.
1
#include "
PhysicsTools/SelectorUtils/interface/CutApplicatorBase.h
"
2
#include "
DataFormats/EgammaCandidates/interface/GsfElectron.h
"
3
4
class
GsfEleMinEcalEtCut
:
public
CutApplicatorBase
{
5
public
:
6
GsfEleMinEcalEtCut
(
const
edm::ParameterSet
&
c
);
7
8
result_type
operator()
(
const
reco::GsfElectronPtr
&)
const
final
;
9
10
double
value
(
const
reco::CandidatePtr
&
cand
)
const
final
;
11
12
CandidateType
candidateType
()
const
final {
return
ELECTRON
; }
13
14
private
:
15
const
double
minEt_
;
16
};
17
18
DEFINE_EDM_PLUGIN
(
CutApplicatorFactory
,
GsfEleMinEcalEtCut
,
"GsfEleMinEcalEtCut"
);
19
20
GsfEleMinEcalEtCut::GsfEleMinEcalEtCut
(
const
edm::ParameterSet
&
c
)
21
:
CutApplicatorBase
(
c
), minEt_(
c
.getParameter<double>(
"minEt"
)) {}
22
23
CutApplicatorBase::result_type
GsfEleMinEcalEtCut::operator()
(
const
reco::GsfElectronPtr
&
cand
)
const
{
24
const
reco::GsfElectron
& ele = *
cand
;
25
const
float
sinTheta = ele.
p
() != 0 ? ele.
pt
() / ele.
p
() : 0.;
26
const
float
et
= ele.
ecalEnergy
() * sinTheta;
27
return
et
>
minEt_
;
28
}
29
30
double
GsfEleMinEcalEtCut::value
(
const
reco::CandidatePtr
&
cand
)
const
{
31
const
reco::GsfElectronPtr
elePtr(
cand
);
32
const
reco::GsfElectron
& ele = *elePtr;
33
const
float
sinTheta = ele.
p
() != 0 ? ele.
pt
() / ele.
p
() : 0.;
34
const
float
et
= ele.
ecalEnergy
() * sinTheta;
35
return
et
;
36
}
reco::LeafCandidate::pt
double pt() const final
transverse momentum
Definition:
LeafCandidate.h:146
watchdog.const
const
Definition:
watchdog.py:83
candidate_functions::CandidateCut::result_type
bool result_type
Definition:
CandidateCut.h:11
reco::GsfElectron
Definition:
GsfElectron.h:35
GsfElectron.h
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition:
PluginFactory.h:124
GsfEleMinEcalEtCut::value
double value(const reco::CandidatePtr &cand) const final
Definition:
GsfEleMinEcalEtCut.cc:30
GsfEleMinEcalEtCut::GsfEleMinEcalEtCut
GsfEleMinEcalEtCut(const edm::ParameterSet &c)
Definition:
GsfEleMinEcalEtCut.cc:20
edm::ParameterSet
Definition:
ParameterSet.h:47
GsfEleMinEcalEtCut::operator()
result_type operator()(const reco::GsfElectronPtr &) const final
Definition:
GsfEleMinEcalEtCut.cc:23
CutApplicatorBase::CandidateType
CandidateType
Definition:
CutApplicatorBase.h:47
edmplugin::PluginFactory
Definition:
PluginFactory.h:34
cand
Definition:
decayParser.h:32
EgHLTOffHistBins_cfi.et
et
Definition:
EgHLTOffHistBins_cfi.py:8
GsfEleMinEcalEtCut::candidateType
CandidateType candidateType() const final
Definition:
GsfEleMinEcalEtCut.cc:12
edm::Ptr< reco::GsfElectron >
CutApplicatorBase
Definition:
CutApplicatorBase.h:45
GsfEleMinEcalEtCut::minEt_
const double minEt_
Definition:
GsfEleMinEcalEtCut.cc:15
reco::LeafCandidate::p
double p() const final
magnitude of momentum vector
Definition:
LeafCandidate.h:123
GsfEleMinEcalEtCut
Definition:
GsfEleMinEcalEtCut.cc:4
CutApplicatorBase.h
c
auto & c
Definition:
CAHitNtupletGeneratorKernelsImpl.h:56
CutApplicatorBase::ELECTRON
Definition:
CutApplicatorBase.h:47
reco::GsfElectron::ecalEnergy
float ecalEnergy() const
Definition:
GsfElectron.h:883
Generated for CMSSW Reference Manual by
1.8.16