Main Page
Namespaces
Classes
Package Documentation
CommonTools
ParticleFlow
interface
ElectronIDPFCandidateSelectorDefinition.h
Go to the documentation of this file.
1
#ifndef CommonTools_ParticleFlow_ElectronIDPFCandidateSelectorDefinition
2
#define CommonTools_ParticleFlow_ElectronIDPFCandidateSelectorDefinition
3
12
#include "
FWCore/Framework/interface/Event.h
"
13
#include "
FWCore/Framework/interface/ConsumesCollector.h
"
14
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
15
#include "
DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h
"
16
#include "
DataFormats/ParticleFlowCandidate/interface/PFCandidate.h
"
17
#include "
DataFormats/Common/interface/ValueMap.h
"
18
#include "
DataFormats/EgammaCandidates/interface/GsfElectron.h
"
19
#include "
CommonTools/ParticleFlow/interface/PFCandidateSelectorDefinition.h
"
20
#include <algorithm>
21
22
namespace
pf2pat
{
23
24
struct
ElectronIDPFCandidateSelectorDefinition
:
public
PFCandidateSelectorDefinition
{
25
ElectronIDPFCandidateSelectorDefinition
(
const
edm::ParameterSet
&
cfg
,
edm::ConsumesCollector
&& iC)
26
:
electronsToken_
(
27
iC.consumes<
reco
::
GsfElectronCollection
>(cfg.getParameter<
edm
::
InputTag
>(
"recoGsfElectrons"
))),
28
electronIdToken_
(iC.consumes<
edm
::ValueMap<
float
> >(cfg.getParameter<
edm
::
InputTag
>(
"electronIdMap"
))) {
29
if
(cfg.
exists
(
"bitsToCheck"
)) {
30
isBitMap_
=
true
;
31
mask_
= 0;
32
if
(cfg.
existsAs
<std::vector<std::string> >(
"bitsToCheck"
)) {
33
std::vector<std::string> strbits = cfg.
getParameter
<std::vector<std::string> >(
"bitsToCheck"
);
34
for
(std::vector<std::string>::const_iterator istrbit = strbits.begin(), estrbit = strbits.end();
35
istrbit != estrbit;
36
++istrbit) {
37
if
(*istrbit ==
"id"
) {
38
mask_
|= 1;
39
}
else
if
(*istrbit ==
"iso"
) {
40
mask_
|= 2;
41
}
else
if
(*istrbit ==
"conv"
) {
42
mask_
|= 4;
43
}
else
if
(*istrbit ==
"ip"
) {
44
mask_
|= 8;
45
}
else
46
throw
cms::Exception
(
"Configuration"
)
47
<<
"ElectronIDPFCandidateSelector: "
48
<<
"bitsToCheck allowed string values are only id(0), iso(1), conv(2), ip(3).\n"
49
<<
"Otherwise, use uint32_t bitmask).\n"
;
50
}
51
}
else
if
(cfg.
existsAs
<uint32_t>(
"bitsToCheck"
)) {
52
mask_
= cfg.
getParameter
<uint32_t>(
"bitsToCheck"
);
53
}
else
{
54
throw
cms::Exception
(
"Configuration"
)
55
<<
"ElectronIDPFCandidateSelector: "
56
<<
"bitsToCheck must be either a vector of strings, or a uint32 bitmask.\n"
;
57
}
58
}
else
{
59
isBitMap_
=
false
;
60
value_
= cfg.
getParameter
<
double
>(
"electronIdCut"
);
61
}
62
}
63
64
void
select
(
const
HandleToCollection
&
hc
,
const
edm::Event
&
e
,
const
edm::EventSetup
&
s
) {
65
selected_
.clear();
66
67
edm::Handle<reco::GsfElectronCollection>
electrons
;
68
e.
getByToken
(
electronsToken_
, electrons);
69
70
edm::Handle<edm::ValueMap<float>
>
electronId
;
71
e.
getByToken
(
electronIdToken_
, electronId);
72
73
unsigned
key
= 0;
74
for
(collection::const_iterator pfc = hc->begin(); pfc != hc->end(); ++pfc, ++
key
) {
75
// Get GsfTrack for matching with reco::GsfElectron objects
76
reco::GsfTrackRef
PfTk = pfc->gsfTrackRef();
77
78
// skip ones without GsfTrack: they won't be matched anyway
79
if
(PfTk.
isNull
())
80
continue
;
81
82
int
match
= -1;
83
// try first the non-ambiguous tracks
84
for
(reco::GsfElectronCollection::const_iterator it = electrons->begin(), ed = electrons->end(); it != ed;
85
++it) {
86
if
(it->gsfTrack() == PfTk) {
87
match = it - electrons->begin();
88
break
;
89
}
90
}
91
// then the ambiguous ones
92
if
(match == -1) {
93
for
(reco::GsfElectronCollection::const_iterator it = electrons->begin(), ed = electrons->end(); it != ed;
94
++it) {
95
if
(
std::count
(it->ambiguousGsfTracksBegin(), it->ambiguousGsfTracksEnd(), PfTk) > 0) {
96
match = it - electrons->begin();
97
break
;
98
}
99
}
100
}
101
// if found, make a GsfElectronRef and read electron id
102
if
(match != -1) {
103
reco::GsfElectronRef
ref(electrons, match);
104
float
eleId = (*electronId)[ref];
105
bool
pass =
false
;
106
if
(
isBitMap_
) {
107
uint32_t thisval = eleId;
108
pass = ((thisval &
mask_
) ==
mask_
);
109
}
else
{
110
pass = (eleId >
value_
);
111
}
112
if
(pass) {
113
selected_
.push_back(
reco::PFCandidate
(*pfc));
114
reco::PFCandidatePtr
ptrToMother(hc, key);
115
selected_
.back().setSourceCandidatePtr(ptrToMother);
116
}
117
}
118
}
119
}
120
121
private
:
122
edm::EDGetTokenT<reco::GsfElectronCollection>
electronsToken_
;
123
edm::EDGetTokenT<edm::ValueMap<float>
>
electronIdToken_
;
124
bool
isBitMap_
;
125
uint32_t
mask_
;
126
double
value_
;
127
};
128
}
// namespace pf2pat
129
130
#endif
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
pf2pat::ElectronIDPFCandidateSelectorDefinition
Selects PFCandidates basing on cuts provided with string cut parser.
Definition:
ElectronIDPFCandidateSelectorDefinition.h:24
KineDebug3::count
void count()
Definition:
KinematicConstrainedVertexUpdatorT.h:21
HLT_2018_cff.InputTag
InputTag
Definition:
HLT_2018_cff.py:78913
edm::ParameterSet::existsAs
bool existsAs(std::string const ¶meterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition:
ParameterSet.h:160
Exception
Definition:
hltDiff.cc:246
pf2pat::ElectronIDPFCandidateSelectorDefinition::ElectronIDPFCandidateSelectorDefinition
ElectronIDPFCandidateSelectorDefinition(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
Definition:
ElectronIDPFCandidateSelectorDefinition.h:25
edm::Event::getByToken
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition:
Event.h:525
edm::Ref< GsfTrackCollection >
Event.h
alignCSCRings.s
s
Definition:
alignCSCRings.py:92
edm::Handle
Definition:
AssociativeIterator.h:49
edm::ParameterSet::exists
bool exists(std::string const ¶meterName) const
checks if a parameter exists
Definition:
ParameterSet.cc:674
MillePedeFileConverter_cfg.e
e
Definition:
MillePedeFileConverter_cfg.py:37
crabWrapper.key
key
Definition:
crabWrapper.py:19
pf2pat::ElectronIDPFCandidateSelectorDefinition::mask_
uint32_t mask_
Definition:
ElectronIDPFCandidateSelectorDefinition.h:125
ValueMap.h
reco::GsfElectronCollection
std::vector< GsfElectron > GsfElectronCollection
collection of GsfElectron objects
Definition:
GsfElectronFwd.h:14
edm::EDGetTokenT< reco::GsfElectronCollection >
topSingleLeptonDQM_PU_cfi.electronId
electronId
when omitted electron plots will be filled w/o cut on electronId
Definition:
topSingleLeptonDQM_PU_cfi.py:39
ParameterSet.h
PFCandidateSelectorDefinition.h
PFCandidate.h
pf2pat::ElectronIDPFCandidateSelectorDefinition::isBitMap_
bool isBitMap_
Definition:
ElectronIDPFCandidateSelectorDefinition.h:124
pf2pat::ElectronIDPFCandidateSelectorDefinition::electronIdToken_
edm::EDGetTokenT< edm::ValueMap< float > > electronIdToken_
Definition:
ElectronIDPFCandidateSelectorDefinition.h:123
edm::EventSetup
Definition:
EventSetup.h:57
edm::Ptr< PFCandidate >
looper.cfg
cfg
Definition:
looper.py:297
edm::Ref::isNull
bool isNull() const
Checks for null.
Definition:
Ref.h:235
GsfElectron.h
pf2pat::ElectronIDPFCandidateSelectorDefinition::select
void select(const HandleToCollection &hc, const edm::Event &e, const edm::EventSetup &s)
Definition:
ElectronIDPFCandidateSelectorDefinition.h:64
pf2pat
Definition:
ElectronIDPFCandidateSelectorDefinition.h:22
pwdgSkimBPark_cfi.electrons
electrons
Definition:
pwdgSkimBPark_cfi.py:6
pf2pat::PFCandidateSelectorDefinition::selected_
container selected_
Definition:
PFCandidateSelectorDefinition.h:33
pf2pat::ElectronIDPFCandidateSelectorDefinition::electronsToken_
edm::EDGetTokenT< reco::GsfElectronCollection > electronsToken_
Definition:
ElectronIDPFCandidateSelectorDefinition.h:122
reco::PFCandidate
Particle reconstructed by the particle flow algorithm.
Definition:
PFCandidate.h:40
reco
fixed size matrix
Definition:
AlignmentAlgorithmBase.h:45
edm
HLT enums.
Definition:
AlignableModifier.h:19
AnalysisDataFormats_SUSYBSMObjects::hc
susybsm::HSCParticleCollection hc
Definition:
classes.h:25
pf2pat::PFCandidateSelectorDefinition
Definition:
PFCandidateSelectorDefinition.h:10
dqmMemoryStats.float
float
Definition:
dqmMemoryStats.py:127
edm::ParameterSet
Definition:
ParameterSet.h:36
match
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition:
Utils.h:10
edm::Event
Definition:
Event.h:72
pf2pat::ElectronIDPFCandidateSelectorDefinition::value_
double value_
Definition:
ElectronIDPFCandidateSelectorDefinition.h:126
PFCandidateFwd.h
ConsumesCollector.h
edm::ConsumesCollector
Definition:
ConsumesCollector.h:39
Generated for CMSSW Reference Manual by
1.8.11