RecoEgamma
ElectronIdentification
plugins
cuts
GsfEleDxyCut.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/VertexReco/interface/VertexFwd.h
"
4
#include "
DataFormats/VertexReco/interface/Vertex.h
"
5
#include "
DataFormats/GsfTrackReco/interface/GsfTrack.h
"
6
7
class
GsfEleDxyCut
:
public
CutApplicatorWithEventContentBase
{
8
public
:
9
GsfEleDxyCut
(
const
edm::ParameterSet
&
c
);
10
11
result_type
operator()
(
const
reco::GsfElectronPtr
&)
const
final
;
12
13
void
setConsumes
(
edm::ConsumesCollector
&)
final
;
14
void
getEventContent
(
const
edm::EventBase
&)
final
;
15
16
double
value
(
const
reco::CandidatePtr
&
cand
)
const
final
;
17
18
CandidateType
candidateType
()
const
final {
return
ELECTRON
; }
19
20
private
:
21
const
double
_dxyCutValueEB
,
_dxyCutValueEE
,
_barrelCutOff
;
22
edm::Handle<reco::VertexCollection>
_vtxs
;
23
};
24
25
DEFINE_EDM_PLUGIN
(
CutApplicatorFactory
,
GsfEleDxyCut
,
"GsfEleDxyCut"
);
26
27
GsfEleDxyCut::GsfEleDxyCut
(
const
edm::ParameterSet
&
c
)
28
:
CutApplicatorWithEventContentBase
(
c
),
29
_dxyCutValueEB(
c
.getParameter<double>(
"dxyCutValueEB"
)),
30
_dxyCutValueEE(
c
.getParameter<double>(
"dxyCutValueEE"
)),
31
_barrelCutOff(
c
.getParameter<double>(
"barrelCutOff"
)) {
32
edm::InputTag
vertextag =
c
.getParameter<
edm::InputTag
>(
"vertexSrc"
);
33
edm::InputTag
vertextagMiniAOD =
c
.getParameter<
edm::InputTag
>(
"vertexSrcMiniAOD"
);
34
contentTags_
.emplace(
"vertices"
, vertextag);
35
contentTags_
.emplace(
"verticesMiniAOD"
, vertextagMiniAOD);
36
}
37
38
void
GsfEleDxyCut::setConsumes
(
edm::ConsumesCollector
&
cc
) {
39
auto
vtcs =
cc
.mayConsume<
reco::VertexCollection
>(
contentTags_
[
"vertices"
]);
40
auto
vtcsMiniAOD =
cc
.mayConsume<
reco::VertexCollection
>(
contentTags_
[
"verticesMiniAOD"
]);
41
contentTokens_
.emplace(
"vertices"
, vtcs);
42
contentTokens_
.emplace(
"verticesMiniAOD"
, vtcsMiniAOD);
43
}
44
45
void
GsfEleDxyCut::getEventContent
(
const
edm::EventBase
&
ev
) {
46
// First try AOD, then go to miniAOD. Use the same Handle since collection class is the same.
47
ev
.getByLabel(
contentTags_
[
"vertices"
],
_vtxs
);
48
if
(!
_vtxs
.
isValid
())
49
ev
.getByLabel(
contentTags_
[
"verticesMiniAOD"
],
_vtxs
);
50
}
51
52
CutApplicatorBase::result_type
GsfEleDxyCut::operator()
(
const
reco::GsfElectronPtr
&
cand
)
const
{
53
const
float
dxyCutValue =
54
(
std::abs
(
cand
->superCluster()->position().eta()) <
_barrelCutOff
?
_dxyCutValueEB
:
_dxyCutValueEE
);
55
56
const
reco::VertexCollection
& vtxs = *
_vtxs
;
57
const
double
dxy
= (!vtxs.empty() ?
cand
->gsfTrack()->dxy(vtxs[0].
position
()) :
cand
->gsfTrack()->dxy());
58
return
std::abs
(
dxy
) < dxyCutValue;
59
}
60
61
double
GsfEleDxyCut::value
(
const
reco::CandidatePtr
&
cand
)
const
{
62
reco::GsfElectronPtr
ele(
cand
);
63
const
reco::VertexCollection
& vtxs = *
_vtxs
;
64
const
double
dxy
= (!vtxs.empty() ? ele->
gsfTrack
()->dxy(vtxs[0].
position
()) : ele->
gsfTrack
()->dxy());
65
return
std::abs
(
dxy
);
66
}
makeMEIFBenchmarkPlots.ev
ev
Definition:
makeMEIFBenchmarkPlots.py:55
CutApplicatorWithEventContentBase
Definition:
CutApplicatorWithEventContentBase.h:19
Vertex.h
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:341
gpuPixelDoublets::cc
uint32_t cc[maxCellsPerHit]
Definition:
gpuFishbone.h:49
edm::Handle< reco::VertexCollection >
reco::VertexCollection
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition:
VertexFwd.h:9
VertexFwd.h
PVValHelper::dxy
Definition:
PVValidationHelpers.h:48
GsfEleDxyCut::_dxyCutValueEB
const double _dxyCutValueEB
Definition:
GsfEleDxyCut.cc:21
HltBtagPostValidation_cff.c
c
Definition:
HltBtagPostValidation_cff.py:31
GsfEleDxyCut::operator()
result_type operator()(const reco::GsfElectronPtr &) const final
Definition:
GsfEleDxyCut.cc:52
reco::GsfElectron::gsfTrack
GsfTrackRef gsfTrack() const override
reference to a GsfTrack
Definition:
GsfElectron.h:156
GsfEleDxyCut::value
double value(const reco::CandidatePtr &cand) const final
Definition:
GsfEleDxyCut.cc:61
GsfEleDxyCut
Definition:
GsfEleDxyCut.cc:7
edmplugin::PluginFactory
Definition:
PluginFactory.h:35
GsfEleDxyCut::_vtxs
edm::Handle< reco::VertexCollection > _vtxs
Definition:
GsfEleDxyCut.cc:22
candidate_functions::CandidateCut::result_type
bool result_type
Definition:
CandidateCut.h:11
GsfEleDxyCut::_dxyCutValueEE
const double _dxyCutValueEE
Definition:
GsfEleDxyCut.cc:21
funct::abs
Abs< T >::type abs(const T &t)
Definition:
Abs.h:22
GsfEleDxyCut::candidateType
CandidateType candidateType() const final
Definition:
GsfEleDxyCut.cc:18
CutApplicatorWithEventContentBase.h
edm::Ptr< reco::GsfElectron >
GsfEleDxyCut::setConsumes
void setConsumes(edm::ConsumesCollector &) final
Definition:
GsfEleDxyCut.cc:38
edm::EventBase
Definition:
EventBase.h:47
GsfElectron.h
GsfTrack.h
GsfEleDxyCut::getEventContent
void getEventContent(const edm::EventBase &) final
Definition:
GsfEleDxyCut.cc:45
GsfEleDxyCut::GsfEleDxyCut
GsfEleDxyCut(const edm::ParameterSet &c)
Definition:
GsfEleDxyCut.cc:27
CutApplicatorWithEventContentBase::contentTokens_
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
Definition:
CutApplicatorWithEventContentBase.h:40
GsfEleDxyCut::_barrelCutOff
const double _barrelCutOff
Definition:
GsfEleDxyCut.cc:21
edm::HandleBase::isValid
bool isValid() const
Definition:
HandleBase.h:70
edm::InputTag
Definition:
InputTag.h:15
position
static int position[264][3]
Definition:
ReadPGInfo.cc:289
CutApplicatorBase::CandidateType
CandidateType
Definition:
CutApplicatorBase.h:47
edm::ParameterSet
Definition:
ParameterSet.h:47
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition:
PluginFactory.h:123
cand
Definition:
decayParser.h:32
edm::ConsumesCollector
Definition:
ConsumesCollector.h:45
Generated for CMSSW Reference Manual by
1.8.14