CMS 3D CMS Logo

ConfigurableTrimmedVertexFinder.cc
Go to the documentation of this file.
3 
4 using namespace reco;
5 
7  const VertexFitter<5> * vf,
8  const VertexUpdator<5> * vu,
10  : theClusterFinder(vf, vu, ve), theVtxFitProbCut(0.01),
11  theTrackCompatibilityToPV(0.05), theTrackCompatibilityToSV(0.01),
12  theMaxNbOfVertices(0)
13 {
14  // default pt cut is 1.5 GeV
15  theFilter.setPtCut(1.5);
16 }
17 
19 {
20  theFilter.setPtCut(s.getParameter<double>("ptCut"));
21  theTrackCompatibilityToPV = s.getParameter<double>("trackCompatibilityToPVcut");
22  theTrackCompatibilityToSV = s.getParameter<double>("trackCompatibilityToSVcut");
23  theVtxFitProbCut = s.getParameter<double>("vtxFitProbCut");
24  theMaxNbOfVertices = s.getParameter<int>("maxNbOfVertices");
25 }
26 
27 
28 std::vector<TransientVertex> ConfigurableTrimmedVertexFinder::vertices(
29  const std::vector<TransientTrack> & tracks) const
30 {
31  std::vector<TransientTrack> remaining;
32 
33  return vertices(tracks, remaining, reco::BeamSpot(), false );
34 
35 }
36 
37 std::vector<TransientVertex> ConfigurableTrimmedVertexFinder::vertices(
38  const std::vector<TransientTrack> & tracks, const reco::BeamSpot & spot ) const
39 {
40  std::vector<TransientTrack> remaining;
41  return vertices ( tracks, remaining, spot, true );
42 }
43 
44 std::vector<TransientVertex> ConfigurableTrimmedVertexFinder::vertices(
45  const std::vector<TransientTrack> & tracks, std::vector<TransientTrack> & unused,
46  const reco::BeamSpot & spot, bool use_spot ) const
47 {
48  resetEvent(tracks);
49  analyseInputTracks(tracks);
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 }
74 
75 
77  const std::vector<TransientTrack> & tracks, std::vector<TransientTrack> & unused,
78  const reco::BeamSpot & spot, bool use_spot ) const
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 }
135 
136 
137 std::vector<TransientVertex>
138 ConfigurableTrimmedVertexFinder::clean(const std::vector<TransientVertex> & candidates) const
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 }
T getParameter(std::string const &) const
std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &tracks) const override
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
float ChiSquaredProbability(double chiSquared, double nrDOF)
void setTrackCompatibilityCut(float cut)
ConfigurableTrimmedVertexFinder(const VertexFitter< 5 > *vf, const VertexUpdator< 5 > *vu, const VertexTrackCompatibilityEstimator< 5 > *ve)
std::vector< TransientVertex > vertexCandidates(const std::vector< reco::TransientTrack > &tracks, std::vector< reco::TransientTrack > &unused, const reco::BeamSpot &spot, bool use_spot) const
void setParameters(const edm::ParameterSet &)
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
fixed size matrix
virtual void analyseVertexCandidates(const std::vector< TransientVertex > &vts) const
void setPtCut(double ptCut)
virtual void resetEvent(const std::vector< reco::TransientTrack > &tracks) const