CMS 3D CMS Logo

ElectronFromPVSelector.cc
Go to the documentation of this file.
1 // Includes
9 
11 
20 
21 #include <algorithm>
22 #include <cmath>
23 #include <memory>
24 #include <vector>
25 
27 // class definition
30 public:
32  void produce(edm::StreamID, edm::Event&, edm::EventSetup const&) const override;
33 
34 private:
35  double max_dxy_;
36  double max_dz_;
39 };
40 
42 // construction
44 
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 
60  iEvent.getByToken(v_recoVertexToken_, vertices);
61 
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 
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