src
RecoEgamma
ElectronIdentification
plugins
cuts
GsfEleConversionVetoCut.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 "
CommonTools/Egamma/interface/ConversionTools.h
"
6
#include "
FWCore/MessageLogger/interface/MessageLogger.h
"
7
8
class
GsfEleConversionVetoCut
:
public
CutApplicatorWithEventContentBase
{
9
public
:
10
GsfEleConversionVetoCut
(
const
edm::ParameterSet
&
c
);
11
12
result_type
operator()
(
const
reco::GsfElectronPtr
&)
const
final
;
13
14
void
setConsumes
(
edm::ConsumesCollector
&)
final
;
15
void
getEventContent
(
const
edm::EventBase
&)
final
;
16
17
double
value
(
const
reco::CandidatePtr
&
cand
)
const
final
;
18
19
CandidateType
candidateType
()
const
final {
return
ELECTRON
; }
20
21
private
:
22
edm::Handle<reco::ConversionCollection>
_convs
;
23
edm::Handle<reco::BeamSpot>
_thebs
;
24
};
25
26
DEFINE_EDM_PLUGIN
(
CutApplicatorFactory
,
GsfEleConversionVetoCut
,
"GsfEleConversionVetoCut"
);
27
28
GsfEleConversionVetoCut::GsfEleConversionVetoCut
(
const
edm::ParameterSet
&
c
) :
CutApplicatorWithEventContentBase
(
c
) {
29
edm::InputTag
conversiontag =
c
.getParameter<
edm::InputTag
>(
"conversionSrc"
);
30
contentTags_
.emplace(
"conversions"
, conversiontag);
31
32
edm::InputTag
conversiontagMiniAOD =
c
.getParameter<
edm::InputTag
>(
"conversionSrcMiniAOD"
);
33
contentTags_
.emplace(
"conversionsMiniAOD"
, conversiontagMiniAOD);
34
35
edm::InputTag
beamspottag =
c
.getParameter<
edm::InputTag
>(
"beamspotSrc"
);
36
contentTags_
.emplace(
"beamspot"
, beamspottag);
37
}
38
39
void
GsfEleConversionVetoCut::setConsumes
(
edm::ConsumesCollector
&
cc
) {
40
auto
convs =
cc
.mayConsume<
reco::ConversionCollection
>(
contentTags_
[
"conversions"
]);
41
auto
convsMiniAOD =
cc
.mayConsume<
reco::ConversionCollection
>(
contentTags_
[
"conversionsMiniAOD"
]);
42
auto
thebs =
cc
.consumes<
reco::BeamSpot
>(
contentTags_
[
"beamspot"
]);
43
contentTokens_
.emplace(
"conversions"
, convs);
44
contentTokens_
.emplace(
"conversionsMiniAOD"
, convsMiniAOD);
45
contentTokens_
.emplace(
"beamspot"
, thebs);
46
}
47
48
void
GsfEleConversionVetoCut::getEventContent
(
const
edm::EventBase
&
ev
) {
49
// First try AOD, then go to miniAOD. Use the same Handle since collection class is the same.
50
ev
.getByLabel(
contentTags_
[
"conversions"
],
_convs
);
51
if
(!
_convs
.
isValid
())
52
ev
.getByLabel(
contentTags_
[
"conversionsMiniAOD"
],
_convs
);
53
54
ev
.getByLabel(
contentTags_
[
"beamspot"
],
_thebs
);
55
}
56
57
CutApplicatorBase::result_type
GsfEleConversionVetoCut::operator()
(
const
reco::GsfElectronPtr
&
cand
)
const
{
58
if
(
_thebs
.
isValid
() &&
_convs
.
isValid
()) {
59
return
!
ConversionTools::hasMatchedConversion
(*
cand
, *
_convs
,
_thebs
->
position
());
60
}
else
{
61
edm::LogWarning
(
"GsfEleConversionVetoCut"
) <<
"Couldn't find a necessary collection, returning true!"
;
62
}
63
return
true
;
64
}
65
66
double
GsfEleConversionVetoCut::value
(
const
reco::CandidatePtr
&
cand
)
const
{
67
reco::GsfElectronPtr
ele(
cand
);
68
if
(
_thebs
.
isValid
() &&
_convs
.
isValid
()) {
69
return
!
ConversionTools::hasMatchedConversion
(*ele, *
_convs
,
_thebs
->
position
());
70
}
else
{
71
edm::LogWarning
(
"GsfEleConversionVetoCut"
) <<
"Couldn't find a necessary collection, returning true!"
;
72
return
true
;
73
}
74
}
makeMEIFBenchmarkPlots.ev
ev
Definition:
makeMEIFBenchmarkPlots.py:55
CutApplicatorWithEventContentBase
Definition:
CutApplicatorWithEventContentBase.h:19
CutApplicatorWithEventContentBase::contentTags_
std::unordered_map< std::string, edm::InputTag > contentTags_
Definition:
CutApplicatorWithEventContentBase.h:35
CutApplicatorBase::ELECTRON
Definition:
CutApplicatorBase.h:47
MessageLogger.h
runTheMatrix.const
const
Definition:
runTheMatrix.py:341
reco::BeamSpot::position
const Point & position() const
position
Definition:
BeamSpot.h:59
GsfEleConversionVetoCut::_thebs
edm::Handle< reco::BeamSpot > _thebs
Definition:
GsfEleConversionVetoCut.cc:23
GsfEleConversionVetoCut::operator()
result_type operator()(const reco::GsfElectronPtr &) const final
Definition:
GsfEleConversionVetoCut.cc:57
gpuPixelDoublets::cc
uint32_t cc[maxCellsPerHit]
Definition:
gpuFishbone.h:49
GsfEleConversionVetoCut::getEventContent
void getEventContent(const edm::EventBase &) final
Definition:
GsfEleConversionVetoCut.cc:48
edm::Handle< reco::ConversionCollection >
Conversion.h
HltBtagPostValidation_cff.c
c
Definition:
HltBtagPostValidation_cff.py:31
reco::ConversionCollection
std::vector< Conversion > ConversionCollection
collectin of Conversion objects
Definition:
ConversionFwd.h:9
GsfEleConversionVetoCut::setConsumes
void setConsumes(edm::ConsumesCollector &) final
Definition:
GsfEleConversionVetoCut.cc:39
GsfEleConversionVetoCut::_convs
edm::Handle< reco::ConversionCollection > _convs
Definition:
GsfEleConversionVetoCut.cc:22
GsfEleConversionVetoCut
Definition:
GsfEleConversionVetoCut.cc:8
edmplugin::PluginFactory
Definition:
PluginFactory.h:35
ConversionTools.h
ConversionTools::hasMatchedConversion
static bool hasMatchedConversion(const reco::GsfElectron &ele, const reco::ConversionCollection &convCol, const math::XYZPoint &beamspot, bool allowCkfMatch=true, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=0)
Definition:
ConversionTools.cc:181
candidate_functions::CandidateCut::result_type
bool result_type
Definition:
CandidateCut.h:11
GsfEleConversionVetoCut::GsfEleConversionVetoCut
GsfEleConversionVetoCut(const edm::ParameterSet &c)
Definition:
GsfEleConversionVetoCut.cc:28
GsfEleConversionVetoCut::candidateType
CandidateType candidateType() const final
Definition:
GsfEleConversionVetoCut.cc:19
CutApplicatorWithEventContentBase.h
edm::Ptr< reco::GsfElectron >
edm::EventBase
Definition:
EventBase.h:47
GsfEleConversionVetoCut::value
double value(const reco::CandidatePtr &cand) const final
Definition:
GsfEleConversionVetoCut.cc:66
GsfElectron.h
CutApplicatorWithEventContentBase::contentTokens_
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
Definition:
CutApplicatorWithEventContentBase.h:40
edm::HandleBase::isValid
bool isValid() const
Definition:
HandleBase.h:70
ConversionFwd.h
edm::InputTag
Definition:
InputTag.h:15
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::LogWarning
Log< level::Warning, false > LogWarning
Definition:
MessageLogger.h:122
reco::BeamSpot
Definition:
BeamSpot.h:21
edm::ConsumesCollector
Definition:
ConsumesCollector.h:45
Generated for CMSSW Reference Manual by
1.8.14