CMS 3D CMS Logo

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

#include <ConfigurableTrimmedVertexFinder.h>

Inheritance diagram for ConfigurableTrimmedVertexFinder:
VertexReconstructor

Public Member Functions

ConfigurableTrimmedVertexFinderclone () const override
 
 ConfigurableTrimmedVertexFinder (const VertexFitter< 5 > *vf, const VertexUpdator< 5 > *vu, const VertexTrackCompatibilityEstimator< 5 > *ve)
 
int maxNbOfVertices () const
 
float ptCut () const
 
void setMaxNbOfVertices (int max)
 
void setParameters (const edm::ParameterSet &)
 
void setPtCut (float cut)
 
void setTrackCompatibilityCut (float cut)
 
void setTrackCompatibilityToSV (float cut)
 
void setVertexFitProbabilityCut (float cut)
 
float trackCompatibilityCut () const
 
float trackCompatibilityToSV () const
 
const TrimmedTrackFiltertrackFilter () const
 
float vertexFitProbabilityCut () const
 
std::vector< TransientVertexvertices (const std::vector< reco::TransientTrack > &tracks) const override
 
std::vector< TransientVertexvertices (const std::vector< reco::TransientTrack > &tracks, const reco::BeamSpot &spot) const override
 
std::vector< TransientVertexvertices (const std::vector< reco::TransientTrack > &tracks, std::vector< reco::TransientTrack > &unused, const reco::BeamSpot &spot, bool use_spot) const
 
 ~ConfigurableTrimmedVertexFinder () override
 
- Public Member Functions inherited from VertexReconstructor
 VertexReconstructor ()
 
virtual std::vector< TransientVertexvertices (const std::vector< reco::TransientTrack > &primaries, const std::vector< reco::TransientTrack > &tracks, const reco::BeamSpot &spot) const
 
virtual ~VertexReconstructor ()
 

Protected Member Functions

virtual void analyseClusterFinder (const std::vector< TransientVertex > &vts, const std::vector< reco::TransientTrack > &remain) const
 
virtual void analyseFoundVertices (const std::vector< TransientVertex > &vts) const
 
virtual void analyseInputTracks (const std::vector< reco::TransientTrack > &tracks) const
 
virtual void analyseVertexCandidates (const std::vector< TransientVertex > &vts) const
 
virtual void resetEvent (const std::vector< reco::TransientTrack > &tracks) const
 

Private Member Functions

std::vector< TransientVertexclean (const std::vector< TransientVertex > &candidates) const
 
std::vector< TransientVertexvertexCandidates (const std::vector< reco::TransientTrack > &tracks, std::vector< reco::TransientTrack > &unused, const reco::BeamSpot &spot, bool use_spot) const
 

Private Attributes

TrimmedVertexFinder theClusterFinder
 
TrimmedTrackFilter theFilter
 
int theMaxNbOfVertices
 
float theTrackCompatibilityToPV
 
float theTrackCompatibilityToSV
 
float theVtxFitProbCut
 

Detailed Description

Algorithm to find a series of distinct vertices among the given set of tracks. The technique is:
1) use TrimmedTrackFilter to select tracks above a certain pT;
2) use TrimmedVertexFinder to split the set of tracks into a cluster of compatible tracks and a set of remaining tracks;
3) repeat 2) with the remaining set, and repeat as long as (a) a cluster of compatible tracks can be found or (b) the maximum number of clusters asked for is reached;
4) reject vertices with a low fit probability.

This algorithm has 5 parameters that can be set at runtime via the corresponding set() methods, or a ParameterSet:

Definition at line 44 of file ConfigurableTrimmedVertexFinder.h.

Constructor & Destructor Documentation

ConfigurableTrimmedVertexFinder::ConfigurableTrimmedVertexFinder ( const VertexFitter< 5 > *  vf,
const VertexUpdator< 5 > *  vu,
const VertexTrackCompatibilityEstimator< 5 > *  ve 
)
ConfigurableTrimmedVertexFinder::~ConfigurableTrimmedVertexFinder ( )
inlineoverride

Definition at line 52 of file ConfigurableTrimmedVertexFinder.h.

References l1t::tracks, and vertices().

52 {}

Member Function Documentation

virtual void ConfigurableTrimmedVertexFinder::analyseClusterFinder ( const std::vector< TransientVertex > &  vts,
const std::vector< reco::TransientTrack > &  remain 
) const
inlineprotectedvirtual

Definition at line 105 of file ConfigurableTrimmedVertexFinder.h.

Referenced by vertexCandidates().

107  {}
virtual void ConfigurableTrimmedVertexFinder::analyseFoundVertices ( const std::vector< TransientVertex > &  vts) const
inlineprotectedvirtual
virtual void ConfigurableTrimmedVertexFinder::analyseInputTracks ( const std::vector< reco::TransientTrack > &  tracks) const
inlineprotectedvirtual

Definition at line 101 of file ConfigurableTrimmedVertexFinder.h.

Referenced by vertices().

102  {}
virtual void ConfigurableTrimmedVertexFinder::analyseVertexCandidates ( const std::vector< TransientVertex > &  vts) const
inlineprotectedvirtual

Definition at line 109 of file ConfigurableTrimmedVertexFinder.h.

Referenced by vertices().

110  {}
std::vector< TransientVertex > ConfigurableTrimmedVertexFinder::clean ( const std::vector< TransientVertex > &  candidates) const
private

Definition at line 138 of file ConfigurableTrimmedVertexFinder.cc.

References ChiSquaredProbability(), mps_fire::i, triggerObjects_cff::sel, and theVtxFitProbCut.

Referenced by analyseFoundVertices(), and vertices().

139 {
140  std::vector<TransientVertex> sel;
141  for (std::vector<TransientVertex>::const_iterator i = candidates.begin();
142  i != candidates.end(); i++) {
143 
144  if (ChiSquaredProbability((*i).totalChiSquared(), (*i).degreesOfFreedom())
145  > theVtxFitProbCut) { sel.push_back(*i); }
146  }
147 
148  return sel;
149 }
float ChiSquaredProbability(double chiSquared, double nrDOF)
ConfigurableTrimmedVertexFinder* ConfigurableTrimmedVertexFinder::clone ( void  ) const
inlineoverridevirtual

Clone method

Implements VertexReconstructor.

Definition at line 93 of file ConfigurableTrimmedVertexFinder.h.

References ConfigurableTrimmedVertexFinder().

Referenced by KalmanTrimmedVertexFinder::KalmanTrimmedVertexFinder().

93  {
94  return new ConfigurableTrimmedVertexFinder(*this);
95  }
ConfigurableTrimmedVertexFinder(const VertexFitter< 5 > *vf, const VertexUpdator< 5 > *vu, const VertexTrackCompatibilityEstimator< 5 > *ve)
int ConfigurableTrimmedVertexFinder::maxNbOfVertices ( ) const
inline
float ConfigurableTrimmedVertexFinder::ptCut ( ) const
inline

Access to parameters

Definition at line 68 of file ConfigurableTrimmedVertexFinder.h.

References TrimmedTrackFilter::ptCut(), and theFilter.

Referenced by KalmanTrimmedVertexFinder::ptCut().

68 { return theFilter.ptCut(); }
double ptCut() const
virtual void ConfigurableTrimmedVertexFinder::resetEvent ( const std::vector< reco::TransientTrack > &  tracks) const
inlineprotectedvirtual

Definition at line 99 of file ConfigurableTrimmedVertexFinder.h.

Referenced by vertices().

99 {}
void ConfigurableTrimmedVertexFinder::setMaxNbOfVertices ( int  max)
inline
void ConfigurableTrimmedVertexFinder::setParameters ( const edm::ParameterSet s)
void ConfigurableTrimmedVertexFinder::setPtCut ( float  cut)
inline
void ConfigurableTrimmedVertexFinder::setTrackCompatibilityCut ( float  cut)
inline
void ConfigurableTrimmedVertexFinder::setTrackCompatibilityToSV ( float  cut)
inline
void ConfigurableTrimmedVertexFinder::setVertexFitProbabilityCut ( float  cut)
inline
float ConfigurableTrimmedVertexFinder::trackCompatibilityCut ( ) const
inline
float ConfigurableTrimmedVertexFinder::trackCompatibilityToSV ( ) const
inline
const TrimmedTrackFilter& ConfigurableTrimmedVertexFinder::trackFilter ( ) const
inline

Definition at line 69 of file ConfigurableTrimmedVertexFinder.h.

References theFilter.

69  {
70  return theFilter;
71  }
std::vector< TransientVertex > ConfigurableTrimmedVertexFinder::vertexCandidates ( const std::vector< reco::TransientTrack > &  tracks,
std::vector< reco::TransientTrack > &  unused,
const reco::BeamSpot spot,
bool  use_spot 
) const
private

Definition at line 76 of file ConfigurableTrimmedVertexFinder.cc.

References analyseClusterFinder(), createfilelist::int, TrimmedVertexFinder::setTrackCompatibilityCut(), theClusterFinder, theMaxNbOfVertices, theTrackCompatibilityToPV, theTrackCompatibilityToSV, l1t::tracks, and TrimmedVertexFinder::vertices().

Referenced by analyseFoundVertices(), and vertices().

79 {
80 
81  std::vector<TransientVertex> cand;
82 
83  std::vector<TransientTrack> remain = tracks;
84 
85  while (true) {
86 
87  float tkCompCut = (cand.empty() ?
90 
91  // std::cout << "PVR:compat cut " << tkCompCut << std::endl;
93  // std::cout << "PVCF:compat cut after setting "
94  // << theClusterFinder.trackCompatibilityCut() << std::endl;
95 
96  std::vector<TransientVertex> newVertices;
97  if ( cand.empty() && use_spot )
98  {
99  newVertices = theClusterFinder.vertices(remain, spot );
100  } else {
101  newVertices = theClusterFinder.vertices(remain);
102  }
103  if (newVertices.empty()) break;
104 
105  analyseClusterFinder(newVertices, remain);
106 
107  for (std::vector<TransientVertex>::const_iterator iv = newVertices.begin();
108  iv != newVertices.end(); iv++) {
109  if ( iv->originalTracks().size() > 1 ) {
110  cand.push_back(*iv);
111  }
112  else {
113  // candidate has too few tracks - get them back into the vector
114  for ( std::vector< TransientTrack >::const_iterator trk
115  = iv->originalTracks().begin();
116  trk != iv->originalTracks().end(); ++trk ) {
117  unused.push_back ( *trk );
118  }
119  }
120  }
121 
122  // when max number of vertices reached, stop
123  if (theMaxNbOfVertices != 0) {
124  if (cand.size() >= (unsigned int) theMaxNbOfVertices) break;
125  }
126  }
127 
128  for (std::vector<TransientTrack>::const_iterator it = remain.begin();
129  it != remain.end(); it++) {
130  unused.push_back(*it);
131  }
132 
133  return cand;
134 }
void setTrackCompatibilityCut(float cut)
std::vector< TransientVertex > vertices(std::vector< reco::TransientTrack > &remain) const
virtual void analyseClusterFinder(const std::vector< TransientVertex > &vts, const std::vector< reco::TransientTrack > &remain) const
float ConfigurableTrimmedVertexFinder::vertexFitProbabilityCut ( ) const
inline
std::vector< TransientVertex > ConfigurableTrimmedVertexFinder::vertices ( const std::vector< reco::TransientTrack > &  ) const
overridevirtual

Reconstruct vertices

Implements VertexReconstructor.

Definition at line 28 of file ConfigurableTrimmedVertexFinder.cc.

Referenced by KalmanTrimmedVertexFinder::vertices(), vertices(), and ~ConfigurableTrimmedVertexFinder().

30 {
31  std::vector<TransientTrack> remaining;
32 
33  return vertices(tracks, remaining, reco::BeamSpot(), false );
34 
35 }
std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &tracks) const override
std::vector< TransientVertex > ConfigurableTrimmedVertexFinder::vertices ( const std::vector< reco::TransientTrack > &  t,
const reco::BeamSpot  
) const
overridevirtual

Reconstruct vertices, exploiting the beamspot constraint for the primary vertex

Reimplemented from VertexReconstructor.

Definition at line 37 of file ConfigurableTrimmedVertexFinder.cc.

References vertices().

39 {
40  std::vector<TransientTrack> remaining;
41  return vertices ( tracks, remaining, spot, true );
42 }
std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &tracks) const override
std::vector< TransientVertex > ConfigurableTrimmedVertexFinder::vertices ( const std::vector< reco::TransientTrack > &  tracks,
std::vector< reco::TransientTrack > &  unused,
const reco::BeamSpot spot,
bool  use_spot 
) const

Definition at line 44 of file ConfigurableTrimmedVertexFinder.cc.

References Vispa.Plugins.EdmBrowser.EdmDataAccessor::all(), analyseFoundVertices(), analyseInputTracks(), analyseVertexCandidates(), clean(), filtered(), resetEvent(), triggerObjects_cff::sel, theFilter, and vertexCandidates().

47 {
50 
51  std::vector<TransientTrack> filtered;
52  for (std::vector<TransientTrack>::const_iterator it = tracks.begin();
53  it != tracks.end(); it++) {
54  if (theFilter(*it)) {
55  filtered.push_back(*it);
56  }
57  else {
58  unused.push_back(*it);
59  }
60  }
61 
62  std::vector<TransientVertex> all = vertexCandidates(filtered, unused,
63  spot, use_spot );
64 
66 
67  std::vector<TransientVertex> sel = clean(all);
68 
70 
71  return sel;
72 
73 }
virtual void analyseFoundVertices(const std::vector< TransientVertex > &vts) const
std::vector< TransientVertex > clean(const std::vector< TransientVertex > &candidates) const
static const TGPicture * filtered(bool iBackgroundIsBlack)
virtual void analyseInputTracks(const std::vector< reco::TransientTrack > &tracks) const
std::vector< TransientVertex > vertexCandidates(const std::vector< reco::TransientTrack > &tracks, std::vector< reco::TransientTrack > &unused, const reco::BeamSpot &spot, bool use_spot) const
virtual void analyseVertexCandidates(const std::vector< TransientVertex > &vts) const
virtual void resetEvent(const std::vector< reco::TransientTrack > &tracks) const

Member Data Documentation

TrimmedVertexFinder ConfigurableTrimmedVertexFinder::theClusterFinder
mutableprivate

Definition at line 135 of file ConfigurableTrimmedVertexFinder.h.

Referenced by vertexCandidates().

TrimmedTrackFilter ConfigurableTrimmedVertexFinder::theFilter
private
int ConfigurableTrimmedVertexFinder::theMaxNbOfVertices
private
float ConfigurableTrimmedVertexFinder::theTrackCompatibilityToPV
private
float ConfigurableTrimmedVertexFinder::theTrackCompatibilityToSV
private
float ConfigurableTrimmedVertexFinder::theVtxFitProbCut
private