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 
26 
28 // class definition
31 public:
32 
34  void produce(edm::StreamID, edm::Event&, edm::EventSetup const&) const override;
35 
36 private:
37  double max_dxy_;
38  double max_dz_;
41 };
42 
43 
45 // construction
47 
49  : max_dxy_{iConfig.getParameter<double>("max_dxy")}
50  , max_dz_{iConfig.getParameter<double>("max_dz")}
51  , v_recoVertexToken_{consumes<std::vector<reco::Vertex>>(iConfig.getParameter<edm::InputTag>("srcVertex"))}
52  , v_recoGsfElectronToken_{consumes<std::vector<reco::GsfElectron>>(iConfig.getParameter<edm::InputTag>("srcElectron"))}
53 {
54  produces<std::vector<reco::GsfElectron>>();
55 }
56 
58 // implementation of member functions
60 
62 {
64  iEvent.getByToken(v_recoVertexToken_, vertices);
65 
67  iEvent.getByToken(v_recoGsfElectronToken_, gsfElectrons);
68 
69  auto goodGsfElectrons = std::make_unique<std::vector<reco::GsfElectron>>();
70 
71  if (!vertices->empty() && !gsfElectrons->empty()) {
72  auto const& pv = vertices->front();
73  std::copy_if(std::cbegin(*gsfElectrons), std::cend(*gsfElectrons), std::back_inserter(*goodGsfElectrons),
74  [this, &pv](auto const& GsfElectron) {
75  return std::abs(GsfElectron.gsfTrack()->dxy(pv.position())) < max_dxy_ &&
76  std::abs(GsfElectron.gsfTrack()->dz(pv.position())) < max_dz_;
77  });
78  }
79 
80  iEvent.put(std::move(goodGsfElectrons));
81 }
82 
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
GsfElectronFromPVSelector(edm::ParameterSet const &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::EDGetTokenT< std::vector< reco::GsfElectron > > v_recoGsfElectronToken_
int iEvent
Definition: GenABIO.cc:230
def pv(vc)
Definition: MetAnalyzer.py:6
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void produce(edm::StreamID, edm::Event &, edm::EventSetup const &) const override
def move(src, dest)
Definition: eostools.py:510
edm::EDGetTokenT< std::vector< reco::Vertex > > v_recoVertexToken_