CMS 3D CMS Logo

VertexMerging.cc
Go to the documentation of this file.
2 
4  : maxFraction(params.getParameter<double>("maxFraction")),
5  minSignificance(params.getParameter<double>("minSignificance")) {}
6 
7 static double computeSharedTracks(const reco::Vertex &pv, const reco::Vertex &sv) {
8  std::set<reco::TrackRef> pvTracks;
9  for (std::vector<reco::TrackBaseRef>::const_iterator iter = pv.tracks_begin(); iter != pv.tracks_end(); iter++) {
10  if (pv.trackWeight(*iter) >= 0.5)
11  pvTracks.insert(iter->castTo<reco::TrackRef>());
12  }
13 
14  unsigned int count = 0, total = 0;
15  for (std::vector<reco::TrackBaseRef>::const_iterator iter = sv.tracks_begin(); iter != sv.tracks_end(); iter++) {
16  if (sv.trackWeight(*iter) >= 0.5) {
17  total++;
18  count += pvTracks.count(iter->castTo<reco::TrackRef>());
19  }
20  }
21 
22  return (double)count / (double)total;
23 }
24 
26  VertexDistance3D dist;
27  reco::VertexCollection recoVertices;
28  for (std::vector<reco::Vertex>::const_iterator sv = secondaryVertices.begin(); sv != secondaryVertices.end(); ++sv) {
29  recoVertices.push_back(*sv);
30  }
31  for (std::vector<reco::Vertex>::iterator sv = recoVertices.begin(); sv != recoVertices.end(); ++sv) {
32  bool shared = false;
33  for (std::vector<reco::Vertex>::iterator sv2 = recoVertices.begin(); sv2 != recoVertices.end(); ++sv2) {
34  double fr = computeSharedTracks(*sv2, *sv);
35  // std::cout << sv2-recoVertices->begin() << " vs " << sv-recoVertices->begin() << " : " << fr << " " << computeSharedTracks(*sv, *sv2) << " sig " << dist.distance(*sv,*sv2).significance() << std::endl;
36  // std::cout << (fr > maxFraction) << " && " << (dist.distance(*sv,*sv2).significance() < 2) << " && " << (sv-sv2!=0) << " && " << (fr >= computeSharedTracks(*sv2, *sv)) << std::endl;
37  if (fr > maxFraction && dist.distance(*sv, *sv2).significance() < minSignificance && sv - sv2 != 0 &&
38  fr >= computeSharedTracks(*sv, *sv2)) {
39  shared = true;
40  // std::cout << "shared " << sv-recoVertices->begin() << " and " << sv2-recoVertices->begin() << " fractions: " << fr << " , " << computeSharedTracks(*sv2, *sv) << " sig: " << dist.distance(*sv,*sv2).significance() << std::endl;
41  }
42  }
43  if (shared) {
44  sv = recoVertices.erase(sv) - 1;
45  }
46  // std::cout << "it = " << sv-recoVertices->begin() << " new size is: " << recoVertices->size() << std::endl;
47  }
48 
49  return recoVertices;
50 }
double minSignificance
Definition: VertexMerging.h:26
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
reco::VertexCollection mergeVertex(reco::VertexCollection &secondaryVertices)
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
double maxFraction
Definition: VertexMerging.h:25
double significance() const
Definition: Measurement1D.h:29
VertexMerging(const edm::ParameterSet &params)
Definition: VertexMerging.cc:3
static double computeSharedTracks(const reco::Vertex &pv, const reco::Vertex &sv)
Definition: VertexMerging.cc:7