Main Page
Namespaces
Classes
Package Documentation
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Pages
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 "
RecoEgamma/EgammaTools/interface/ConversionTools.h
"
6
7
8
class
GsfEleConversionVetoCut
:
public
CutApplicatorWithEventContentBase
{
9
public
:
10
GsfEleConversionVetoCut
(
const
edm::ParameterSet
&
c
);
11
12
result_type
operator()
(
const
reco::GsfElectronPtr
&)
const
override
final
;
13
14
void
setConsumes
(
edm::ConsumesCollector
&)
override
final
;
15
void
getEventContent
(
const
edm::EventBase
&)
override
final
;
16
17
CandidateType
candidateType
()
const
override final {
18
return
ELECTRON
;
19
}
20
21
private
:
22
edm::Handle<reco::ConversionCollection>
_convs
;
23
edm::Handle<reco::BeamSpot>
_thebs
;
24
};
25
26
DEFINE_EDM_PLUGIN
(
CutApplicatorFactory
,
27
GsfEleConversionVetoCut
,
28
"GsfEleConversionVetoCut"
);
29
30
GsfEleConversionVetoCut::GsfEleConversionVetoCut
(
const
edm::ParameterSet
&
c
) :
31
CutApplicatorWithEventContentBase
(c) {
32
edm::InputTag
conversiontag = c.
getParameter
<
edm::InputTag
>(
"conversionSrc"
);
33
contentTags_
.emplace(
"conversions"
,conversiontag);
34
edm::InputTag
beamspottag = c.
getParameter
<
edm::InputTag
>(
"beamspotSrc"
);
35
contentTags_
.emplace(
"beamspot"
,beamspottag);
36
}
37
38
void
GsfEleConversionVetoCut::setConsumes
(
edm::ConsumesCollector
& cc) {
39
auto
convs =
40
cc.
consumes
<
reco::ConversionCollection
>(
contentTags_
[
"conversions"
]);
41
auto
thebs = cc.
consumes
<
reco::BeamSpot
>(
contentTags_
[
"beamspot"
]);
42
contentTokens_
.emplace(
"conversions"
,convs);
43
contentTokens_
.emplace(
"beamspot"
,thebs);
44
}
45
46
void
GsfEleConversionVetoCut::getEventContent
(
const
edm::EventBase
& ev) {
47
ev.
getByLabel
(
contentTags_
[
"conversions"
],
_convs
);
48
ev.
getByLabel
(
contentTags_
[
"beamspot"
],
_thebs
);
49
}
50
51
CutApplicatorBase::result_type
52
GsfEleConversionVetoCut::
53
operator()
(
const
reco::GsfElectronPtr
& cand)
const
{
54
if
(
_thebs
.
isValid
() &&
_convs
.
isValid
() ) {
55
return
!
ConversionTools::hasMatchedConversion
(*cand,
_convs
,
56
_thebs
->position());
57
}
else
{
58
edm::LogWarning
(
"GsfEleConversionVetoCut"
)
59
<<
"Couldn't find a necessary collection, returning true!"
;
60
}
61
return
true
;
62
}
edm::ConsumesCollector::consumes
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
Definition:
ConsumesCollector.h:41
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
CutApplicatorWithEventContentBase
Definition:
CutApplicatorWithEventContentBase.h:18
CutApplicatorBase::ELECTRON
Definition:
CutApplicatorBase.h:48
GsfEleConversionVetoCut::_thebs
edm::Handle< reco::BeamSpot > _thebs
Definition:
GsfEleConversionVetoCut.cc:23
GsfEleConversionVetoCut::candidateType
CandidateType candidateType() const overridefinal
Definition:
GsfEleConversionVetoCut.cc:17
GsfEleConversionVetoCut::setConsumes
void setConsumes(edm::ConsumesCollector &) overridefinal
Definition:
GsfEleConversionVetoCut.cc:38
CutApplicatorWithEventContentBase::contentTokens_
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
Definition:
CutApplicatorWithEventContentBase.h:40
edm::Handle< reco::ConversionCollection >
edm::LogWarning
Definition:
MessageLogger.h:140
Conversion.h
CutApplicatorWithEventContentBase::contentTags_
std::unordered_map< std::string, edm::InputTag > contentTags_
Definition:
CutApplicatorWithEventContentBase.h:35
reco::ConversionCollection
std::vector< Conversion > ConversionCollection
collectin of Conversion objects
Definition:
ConversionFwd.h:9
GsfEleConversionVetoCut::_convs
edm::Handle< reco::ConversionCollection > _convs
Definition:
GsfEleConversionVetoCut.cc:22
GsfEleConversionVetoCut
Definition:
GsfEleConversionVetoCut.cc:8
edmplugin::PluginFactory
Definition:
PluginFactory.h:31
ConversionTools.h
GsfEleConversionVetoCut::GsfEleConversionVetoCut
GsfEleConversionVetoCut(const edm::ParameterSet &c)
Definition:
GsfEleConversionVetoCut.cc:30
CutApplicatorWithEventContentBase.h
edm::Ptr< reco::GsfElectron >
edm::HandleBase::isValid
bool isValid() const
Definition:
HandleBase.h:76
ConversionTools::hasMatchedConversion
static bool hasMatchedConversion(const reco::GsfElectron &ele, const edm::Handle< 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:149
edm::EventBase
Definition:
EventBase.h:45
GsfElectron.h
GsfEleConversionVetoCut::getEventContent
void getEventContent(const edm::EventBase &) overridefinal
Definition:
GsfEleConversionVetoCut.cc:46
GsfEleConversionVetoCut::operator()
result_type operator()(const reco::GsfElectronPtr &) const overridefinal
Definition:
GsfEleConversionVetoCut.cc:53
compareJSON.const
string const
Definition:
compareJSON.py:14
trackerHits.c
tuple c
Definition:
trackerHits.py:26
ConversionFwd.h
edm::EventBase::getByLabel
bool getByLabel(InputTag const &, Handle< T > &) const
Definition:
EventBase.h:87
edm::InputTag
Definition:
InputTag.h:17
CutApplicatorBase::CandidateType
CandidateType
Definition:
CutApplicatorBase.h:47
edm::ParameterSet
Definition:
ParameterSet.h:35
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition:
PluginFactory.h:101
reco::BeamSpot
Definition:
BeamSpot.h:22
edm::ConsumesCollector
Definition:
ConsumesCollector.h:32
Generated for CMSSW Reference Manual by
1.8.5