CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
VertexMerging Class Reference

#include <VertexMerging.h>

Public Member Functions

reco::VertexCollection mergeVertex (reco::VertexCollection &secondaryVertices)
 
 VertexMerging (const edm::ParameterSet &params)
 

Private Member Functions

bool trackFilter (const reco::TrackRef &track) const
 

Private Attributes

double maxFraction
 
double minSignificance
 

Detailed Description

Definition at line 17 of file VertexMerging.h.

Constructor & Destructor Documentation

VertexMerging::VertexMerging ( const edm::ParameterSet params)

Definition at line 3 of file VertexMerging.cc.

4  : maxFraction(params.getParameter<double>("maxFraction")),
5  minSignificance(params.getParameter<double>("minSignificance")) {}
T getParameter(std::string const &) const
double minSignificance
Definition: VertexMerging.h:27
double maxFraction
Definition: VertexMerging.h:26

Member Function Documentation

reco::VertexCollection VertexMerging::mergeVertex ( reco::VertexCollection secondaryVertices)

Definition at line 25 of file VertexMerging.cc.

References computeSharedTracks(), VertexDistance3D::distance(), maxFraction, minSignificance, Measurement1D::significance(), and pfDeepBoostedJetPreprocessParams_cfi::sv.

25  {
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:27
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
double significance() const
Definition: Measurement1D.h:29
double maxFraction
Definition: VertexMerging.h:26
static double computeSharedTracks(const reco::Vertex &pv, const reco::Vertex &sv)
Definition: VertexMerging.cc:7
bool VertexMerging::trackFilter ( const reco::TrackRef track) const
private

Member Data Documentation

double VertexMerging::maxFraction
private

Definition at line 26 of file VertexMerging.h.

Referenced by mergeVertex().

double VertexMerging::minSignificance
private

Definition at line 27 of file VertexMerging.h.

Referenced by mergeVertex().