Main Page
Namespaces
Classes
Package Documentation
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Pages
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/ConsumesCollector.h
"
13
#include "
DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h
"
14
#include "
DataFormats/ParticleFlowCandidate/interface/PFCandidate.h
"
15
#include "
DataFormats/Common/interface/ValueMap.h
"
16
#include "
DataFormats/EgammaCandidates/interface/GsfElectron.h
"
17
#include "
CommonTools/ParticleFlow/interface/PFCandidateSelectorDefinition.h
"
18
#include <algorithm>
19
20
namespace
pf2pat {
21
22
struct
ElectronIDPFCandidateSelectorDefinition
:
public
PFCandidateSelectorDefinition
{
23
24
ElectronIDPFCandidateSelectorDefinition
(
const
edm::ParameterSet
&
cfg
,
edm::ConsumesCollector
&& iC ) :
25
electronsToken_
( iC.consumes<
reco
::
GsfElectronCollection
>( cfg.getParameter< edm::
InputTag
>(
"recoGsfElectrons"
) ) ),
26
electronIdToken_
( iC.consumes<edm::ValueMap<float> >( cfg.getParameter< edm::
InputTag
>(
"electronIdMap"
) ) )
27
{
28
if
(cfg.
exists
(
"bitsToCheck"
)) {
29
isBitMap_
=
true
;
30
mask_
= 0;
31
if
(cfg.
existsAs
<std::vector<std::string> >(
"bitsToCheck"
)) {
32
std::vector<std::string> strbits = cfg.
getParameter
<std::vector<std::string> >(
"bitsToCheck"
);
33
for
(std::vector<std::string>::const_iterator istrbit = strbits.begin(), estrbit = strbits.end();
34
istrbit != estrbit; ++istrbit) {
35
if
(*istrbit ==
"id"
) {
mask_
|= 1; }
36
else
if
(*istrbit ==
"iso"
) {
mask_
|= 2; }
37
else
if
(*istrbit ==
"conv"
) {
mask_
|= 4; }
38
else
if
(*istrbit ==
"ip"
) {
mask_
|= 8; }
39
else
throw
cms::Exception
(
"Configuration"
) <<
"ElectronIDPFCandidateSelector: "
<<
40
"bitsToCheck allowed string values are only id(0), iso(1), conv(2), ip(3).\n"
<<
41
"Otherwise, use uint32_t bitmask).\n"
;
42
}
43
}
else
if
(cfg.
existsAs
<uint32_t>(
"bitsToCheck"
)) {
44
mask_
= cfg.
getParameter
<uint32_t>(
"bitsToCheck"
);
45
}
else
{
46
throw
cms::Exception
(
"Configuration"
) <<
"ElectronIDPFCandidateSelector: "
<<
47
"bitsToCheck must be either a vector of strings, or a uint32 bitmask.\n"
;
48
}
49
}
else
{
50
isBitMap_
=
false
;
51
value_
= cfg.
getParameter
<
double
>(
"electronIdCut"
);
52
}
53
}
54
55
void
select
(
const
HandleToCollection
&
hc
,
56
const
edm::Event
&
e
,
57
const
edm::EventSetup
&
s
) {
58
selected_
.clear();
59
60
edm::Handle<reco::GsfElectronCollection>
electrons
;
61
e.
getByToken
(
electronsToken_
, electrons);
62
63
edm::Handle<edm::ValueMap<float>
> electronId;
64
e.
getByToken
(
electronIdToken_
, electronId);
65
66
unsigned
key
=0;
67
for
( collection::const_iterator pfc = hc->begin();
68
pfc != hc->end(); ++pfc, ++
key
) {
69
70
// Get GsfTrack for matching with reco::GsfElectron objects
71
reco::GsfTrackRef
PfTk = pfc->gsfTrackRef();
72
73
// skip ones without GsfTrack: they won't be matched anyway
74
if
(PfTk.
isNull
())
continue
;
75
76
int
match
= -1;
77
// try first the non-ambiguous tracks
78
for
(reco::GsfElectronCollection::const_iterator it = electrons->begin(), ed = electrons->end(); it != ed; ++it) {
79
if
(it->gsfTrack() == PfTk) { match = it - electrons->begin();
break
; }
80
}
81
// then the ambiguous ones
82
if
(match == -1) {
83
for
(reco::GsfElectronCollection::const_iterator it = electrons->begin(), ed = electrons->end(); it != ed; ++it) {
84
if
(
std::count
(it->ambiguousGsfTracksBegin(), it->ambiguousGsfTracksEnd(), PfTk) > 0) {
85
match = it - electrons->begin();
break
;
86
}
87
}
88
}
89
// if found, make a GsfElectronRef and read electron id
90
if
(match != -1) {
91
reco::GsfElectronRef
ref(electrons,match);
92
float
eleId = (*electronId)[ref];
93
bool
pass =
false
;
94
if
(
isBitMap_
) {
95
uint32_t thisval = eleId;
96
pass = ((thisval &
mask_
) ==
mask_
);
97
}
else
{
98
pass = (eleId >
value_
);
99
}
100
if
(pass) {
101
selected_
.push_back(
reco::PFCandidate
(*pfc) );
102
reco::PFCandidatePtr
ptrToMother( hc, key );
103
selected_
.back().setSourceCandidatePtr( ptrToMother );
104
}
105
}
106
}
107
}
108
109
private
:
110
edm::EDGetTokenT<reco::GsfElectronCollection>
electronsToken_
;
111
edm::EDGetTokenT<edm::ValueMap<float>
>
electronIdToken_
;
112
bool
isBitMap_
;
113
uint32_t
mask_
;
114
double
value_
;
115
};
116
}
117
118
#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:22
KineDebug3::count
void count()
Definition:
KinematicConstrainedVertexUpdatorT.h:20
looper.cfg
tuple cfg
Definition:
looper.py:293
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:186
Exception
Definition:
hltDiff.cc:323
pf2pat::ElectronIDPFCandidateSelectorDefinition::ElectronIDPFCandidateSelectorDefinition
ElectronIDPFCandidateSelectorDefinition(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
Definition:
ElectronIDPFCandidateSelectorDefinition.h:24
edm::Event::getByToken
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition:
Event.h:462
edm::Ref< GsfTrackCollection >
edm::Handle
Definition:
AssociativeIterator.h:47
edm::ParameterSet::exists
bool exists(std::string const ¶meterName) const
checks if a parameter exists
Definition:
ParameterSet.cc:760
pf2pat::ElectronIDPFCandidateSelectorDefinition::mask_
uint32_t mask_
Definition:
ElectronIDPFCandidateSelectorDefinition.h:113
dt_dqm_sourceclient_common_cff.reco
tuple reco
Definition:
dt_dqm_sourceclient_common_cff.py:107
ValueMap.h
HI_PhotonSkim_cff.electrons
tuple electrons
Definition:
HI_PhotonSkim_cff.py:77
reco::GsfElectronCollection
std::vector< GsfElectron > GsfElectronCollection
collection of GsfElectron objects
Definition:
GsfElectronFwd.h:14
edm::EDGetTokenT< reco::GsfElectronCollection >
PFCandidateSelectorDefinition.h
PFCandidate.h
pf2pat::ElectronIDPFCandidateSelectorDefinition::isBitMap_
bool isBitMap_
Definition:
ElectronIDPFCandidateSelectorDefinition.h:112
pf2pat::ElectronIDPFCandidateSelectorDefinition::electronIdToken_
edm::EDGetTokenT< edm::ValueMap< float > > electronIdToken_
Definition:
ElectronIDPFCandidateSelectorDefinition.h:111
edm::EventSetup
Definition:
EventSetup.h:45
edm::Ptr< PFCandidate >
relval_steps.key
string key
FastSim: produces sample of signal events, overlayed with premixed minbias events.
Definition:
relval_steps.py:631
edm::Ref::isNull
bool isNull() const
Checks for null.
Definition:
Ref.h:249
GsfElectron.h
pf2pat::ElectronIDPFCandidateSelectorDefinition::select
void select(const HandleToCollection &hc, const edm::Event &e, const edm::EventSetup &s)
Definition:
ElectronIDPFCandidateSelectorDefinition.h:55
alignCSCRings.s
list s
Definition:
alignCSCRings.py:91
pf2pat::PFCandidateSelectorDefinition::selected_
container selected_
Definition:
PFCandidateSelectorDefinition.h:35
alignCSCRings.e
list e
Definition:
alignCSCRings.py:90
HLT_FULL_cff.InputTag
tuple InputTag
Definition:
HLT_FULL_cff.py:69759
pf2pat::ElectronIDPFCandidateSelectorDefinition::electronsToken_
edm::EDGetTokenT< reco::GsfElectronCollection > electronsToken_
Definition:
ElectronIDPFCandidateSelectorDefinition.h:110
reco::PFCandidate
Particle reconstructed by the particle flow algorithm.
Definition:
PFCandidate.h:39
AnalysisDataFormats_SUSYBSMObjects::hc
susybsm::HSCParticleCollection hc
Definition:
classes.h:25
pf2pat::PFCandidateSelectorDefinition
Definition:
PFCandidateSelectorDefinition.h:11
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:65
pf2pat::ElectronIDPFCandidateSelectorDefinition::value_
double value_
Definition:
ElectronIDPFCandidateSelectorDefinition.h:114
PFCandidateFwd.h
ConsumesCollector.h
edm::ConsumesCollector
Definition:
ConsumesCollector.h:39
Generated for CMSSW Reference Manual by
1.8.5