Validation
RecoTau
plugins
ElectronFromPVSelector.cc
Go to the documentation of this file.
1
// Includes
4
#include "
FWCore/Framework/interface/Event.h
"
5
#include "
FWCore/Framework/interface/MakerMacros.h
"
6
#include "
FWCore/Framework/interface/global/EDProducer.h
"
7
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
8
#include "
FWCore/Utilities/interface/InputTag.h
"
9
10
#include "
CommonTools/Utils/interface/StringCutObjectSelector.h
"
11
12
#include "
DataFormats/EgammaCandidates/interface/GsfElectron.h
"
13
#include "
DataFormats/EgammaCandidates/interface/GsfElectronFwd.h
"
14
#include "
DataFormats/GsfTrackReco/interface/GsfTrack.h
"
15
#include "
DataFormats/GsfTrackReco/interface/GsfTrackExtra.h
"
16
#include "
DataFormats/GsfTrackReco/interface/GsfTrackExtraFwd.h
"
17
#include "
DataFormats/GsfTrackReco/interface/GsfTrackFwd.h
"
18
#include "
DataFormats/VertexReco/interface/Vertex.h
"
19
#include "
DataFormats/VertexReco/interface/VertexFwd.h
"
20
21
#include <algorithm>
22
#include <cmath>
23
#include <memory>
24
#include <vector>
25
27
// class definition
29
class
GsfElectronFromPVSelector
:
public
edm::global::EDProducer
<> {
30
public
:
31
explicit
GsfElectronFromPVSelector
(
edm::ParameterSet
const
&);
32
void
produce
(
edm::StreamID
,
edm::Event
&,
edm::EventSetup
const
&)
const override
;
33
34
private
:
35
double
max_dxy_
;
36
double
max_dz_
;
37
edm::EDGetTokenT<std::vector<reco::Vertex>
>
v_recoVertexToken_
;
38
edm::EDGetTokenT<std::vector<reco::GsfElectron>
>
v_recoGsfElectronToken_
;
39
};
40
42
// construction
44
45
GsfElectronFromPVSelector::GsfElectronFromPVSelector
(
edm::ParameterSet
const
& iConfig)
46
: max_dxy_{iConfig.
getParameter
<
double
>(
"max_dxy"
)},
47
max_dz_{iConfig.getParameter<
double
>(
"max_dz"
)},
48
v_recoVertexToken_{consumes<std::vector<reco::Vertex>>(iConfig.getParameter<
edm::InputTag
>(
"srcVertex"
))},
49
v_recoGsfElectronToken_{
50
consumes<std::vector<reco::GsfElectron>>(iConfig.getParameter<
edm::InputTag
>(
"srcElectron"
))} {
51
produces<std::vector<reco::GsfElectron>>();
52
}
53
55
// implementation of member functions
57
58
void
GsfElectronFromPVSelector::produce
(
edm::StreamID
,
edm::Event
&
iEvent
,
edm::EventSetup
const
&)
const
{
59
edm::Handle<std::vector<reco::Vertex>
>
vertices
;
60
iEvent
.getByToken(
v_recoVertexToken_
,
vertices
);
61
62
edm::Handle<std::vector<reco::GsfElectron>
>
gsfElectrons
;
63
iEvent
.getByToken(
v_recoGsfElectronToken_
,
gsfElectrons
);
64
65
auto
goodGsfElectrons = std::make_unique<std::vector<reco::GsfElectron>>();
66
67
if
(!
vertices
->empty() && !
gsfElectrons
->empty()) {
68
auto
const
&
pv
=
vertices
->front();
69
std::copy_if(std::cbegin(*
gsfElectrons
),
70
std::cend(*
gsfElectrons
),
71
std::back_inserter(*goodGsfElectrons),
72
[
this
, &
pv
](
auto
const
& GsfElectron) {
73
return
std::abs
(GsfElectron.gsfTrack()->dxy(
pv
.position())) <
max_dxy_
&&
74
std::abs
(GsfElectron.gsfTrack()->dz(
pv
.position())) <
max_dz_
;
75
});
76
}
77
78
iEvent
.put(
std::move
(goodGsfElectrons));
79
}
80
81
DEFINE_FWK_MODULE
(
GsfElectronFromPVSelector
);
GsfElectronFromPVSelector::GsfElectronFromPVSelector
GsfElectronFromPVSelector(edm::ParameterSet const &)
Definition:
ElectronFromPVSelector.cc:45
edm::StreamID
Definition:
StreamID.h:30
GsfTrackExtra.h
edm::EDGetTokenT
Definition:
EDGetToken.h:33
GsfElectronFromPVSelector::produce
void produce(edm::StreamID, edm::Event &, edm::EventSetup const &) const override
Definition:
ElectronFromPVSelector.cc:58
electronIsolatorFromEffectiveArea_cfi.gsfElectrons
gsfElectrons
Definition:
electronIsolatorFromEffectiveArea_cfi.py:4
edm::Handle
Definition:
AssociativeIterator.h:50
GsfElectronFromPVSelector::v_recoVertexToken_
edm::EDGetTokenT< std::vector< reco::Vertex > > v_recoVertexToken_
Definition:
ElectronFromPVSelector.cc:37
GsfElectronFromPVSelector::v_recoGsfElectronToken_
edm::EDGetTokenT< std::vector< reco::GsfElectron > > v_recoGsfElectronToken_
Definition:
ElectronFromPVSelector.cc:38
MakerMacros.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition:
MakerMacros.h:16
GsfElectronFromPVSelector::max_dz_
double max_dz_
Definition:
ElectronFromPVSelector.cc:36
GsfElectron.h
edm::global::EDProducer
Definition:
EDProducer.h:32
Vertex.h
GsfElectronFwd.h
edm::ParameterSet
Definition:
ParameterSet.h:36
Event.h
GsfElectronFromPVSelector::max_dxy_
double max_dxy_
Definition:
ElectronFromPVSelector.cc:35
GsfTrackExtraFwd.h
MetAnalyzer.pv
def pv(vc)
Definition:
MetAnalyzer.py:7
iEvent
int iEvent
Definition:
GenABIO.cc:224
GsfTrack.h
edm::EventSetup
Definition:
EventSetup.h:57
InputTag.h
VertexFwd.h
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
eostools.move
def move(src, dest)
Definition:
eostools.py:511
StringCutObjectSelector.h
GsfTrackFwd.h
funct::abs
Abs< T >::type abs(const T &t)
Definition:
Abs.h:22
ParameterSet.h
EDProducer.h
edm::Event
Definition:
Event.h:73
edm::InputTag
Definition:
InputTag.h:15
GsfElectronFromPVSelector
Definition:
ElectronFromPVSelector.cc:29
pwdgSkimBPark_cfi.vertices
vertices
Definition:
pwdgSkimBPark_cfi.py:7
Generated for CMSSW Reference Manual by
1.8.16