CMS 3D CMS Logo

TrackFromPVSelector.cc
Go to the documentation of this file.
1 // Includes
9 
11 
16 
17 #include <algorithm>
18 #include <numeric>
19 #include <memory>
20 #include <vector>
21 
22 
24 // class definition
27 public:
28 
29  explicit TrackFromPVSelector(edm::ParameterSet const& iConfig);
30 
31  void produce(edm::StreamID, edm::Event& iEvent, edm::EventSetup const& iSetup) const override;
32 
33 private:
34 
35  double max_dxy_;
36  double max_dz_;
39 };
40 
41 
43 // construction/destruction
45 
46 //______________________________________________________________________________
48  : max_dxy_{iConfig.getParameter<double>("max_dxy")}
49  , max_dz_{iConfig.getParameter<double>("max_dz")}
50  , v_recoVertexToken_{consumes<std::vector<reco::Vertex>>(iConfig.getParameter<edm::InputTag>("srcVertex"))}
51  , v_recoTrackToken_{consumes<std::vector<reco::Track>>(iConfig.getParameter<edm::InputTag>("srcTrack"))}
52 {
53  produces<std::vector<reco::Track>>();
54 }
55 
57 // implementation of member functions
59 
60 //______________________________________________________________________________
62 {
64  iEvent.getByToken(v_recoVertexToken_, vertices);
65 
67  iEvent.getByToken(v_recoTrackToken_, tracks);
68 
69  auto goodTracks = std::make_unique<std::vector<reco::Track>>();
70  if (!vertices->empty() && !tracks->empty()) {
71  auto const& vtxPos = vertices->front().position();
72  std::copy_if(std::cbegin(*tracks), std::cend(*tracks), std::back_inserter(*goodTracks),
73  [this, &vtxPos](auto const& track) {
74  return std::abs(track.dxy(vtxPos)) < max_dxy_ && std::abs(track.dz(vtxPos)) < max_dz_;
75  });
76  }
77  iEvent.put(std::move(goodTracks));
78 }
79 
T getParameter(std::string const &) const
edm::EDGetTokenT< std::vector< reco::Vertex > > v_recoVertexToken_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
edm::EDGetTokenT< std::vector< reco::Track > > v_recoTrackToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void produce(edm::StreamID, edm::Event &iEvent, edm::EventSetup const &iSetup) const override
TrackFromPVSelector(edm::ParameterSet const &iConfig)
int iEvent
Definition: GenABIO.cc:230
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
def move(src, dest)
Definition: eostools.py:510