CMS 3D CMS Logo

DoubleVertexFilter.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <set>
3 
9 
15 
17  public:
19 
20  void produce(edm::StreamID, edm::Event &event, const edm::EventSetup &es) const override;
21 
22  private:
23  bool trackFilter(const reco::TrackRef &track) const;
24 
27  double maxFraction;
28 };
29 
31  maxFraction(params.getParameter<double>("maxFraction"))
32 {
33  token_primaryVertex = consumes<reco::VertexCollection>(params.getParameter<edm::InputTag>("primaryVertices"));
34  token_secondaryVertex = consumes<reco::VertexCollection>(params.getParameter<edm::InputTag>("secondaryVertices"));
35  produces<reco::VertexCollection>();
36 }
37 
38 static double computeSharedTracks(const reco::Vertex &pv,
39  const reco::Vertex &sv)
40 {
41  std::set<reco::TrackRef> pvTracks;
42  for(std::vector<reco::TrackBaseRef>::const_iterator iter = pv.tracks_begin();
43  iter != pv.tracks_end(); iter++) {
44  if (pv.trackWeight(*iter) >= 0.5)
45  pvTracks.insert(iter->castTo<reco::TrackRef>());
46  }
47 
48  unsigned int count = 0, total = 0;
49  for(std::vector<reco::TrackBaseRef>::const_iterator iter = sv.tracks_begin();
50  iter != sv.tracks_end(); iter++) {
51  if (sv.trackWeight(*iter) >= 0.5) {
52  total++;
53  count += pvTracks.count(iter->castTo<reco::TrackRef>());
54  }
55  }
56 
57  return (double)count / (double)total;
58 }
59 
61 {
62  using namespace reco;
63 
65  event.getByToken(token_primaryVertex, primaryVertices);
66 
68  event.getByToken(token_secondaryVertex, secondaryVertices);
69 
70  std::vector<reco::Vertex>::const_iterator pv = primaryVertices->begin();
71 
72  auto recoVertices = std::make_unique<VertexCollection>();
73  for(std::vector<reco::Vertex>::const_iterator sv = secondaryVertices->begin();
74  sv != secondaryVertices->end(); ++sv) {
75  if (computeSharedTracks(*pv, *sv) > maxFraction)
76  continue;
77 
78  recoVertices->push_back(*sv);
79  }
80 
81  event.put(std::move(recoVertices));
82 }
83 
T getParameter(std::string const &) const
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:81
bool trackFilter(const reco::TrackRef &track) const
edm::EDGetTokenT< reco::VertexCollection > token_primaryVertex
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static double computeSharedTracks(const reco::Vertex &pv, const reco::Vertex &sv)
def pv(vc)
Definition: MetAnalyzer.py:7
float trackWeight(const TREF &r) const
returns the weight with which a Track has contributed to the vertex-fit.
Definition: Vertex.h:81
void produce(edm::StreamID, edm::Event &event, const edm::EventSetup &es) const override
DoubleVertexFilter(const edm::ParameterSet &params)
edm::EDGetTokenT< reco::VertexCollection > token_secondaryVertex
primaryVertices
Definition: jets_cff.py:24
fixed size matrix
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:76
def move(src, dest)
Definition: eostools.py:511
Definition: event.py:1