CMS 3D CMS Logo

ConfigurableTrimmedVertexFinder.cc
Go to the documentation of this file.
3 
4 using namespace reco;
5 
7  const VertexUpdator<5>* vu,
9  : theClusterFinder(vf, vu, ve),
10  theVtxFitProbCut(0.01),
11  theTrackCompatibilityToPV(0.05),
12  theTrackCompatibilityToSV(0.01),
13  theMaxNbOfVertices(0) {
14  // default pt cut is 1.5 GeV
15  theFilter.setPtCut(1.5);
16 }
17 
19  theFilter.setPtCut(s.getParameter<double>("ptCut"));
20  theTrackCompatibilityToPV = s.getParameter<double>("trackCompatibilityToPVcut");
21  theTrackCompatibilityToSV = s.getParameter<double>("trackCompatibilityToSVcut");
22  theVtxFitProbCut = s.getParameter<double>("vtxFitProbCut");
23  theMaxNbOfVertices = s.getParameter<int>("maxNbOfVertices");
24 }
25 
26 std::vector<TransientVertex> ConfigurableTrimmedVertexFinder::vertices(const std::vector<TransientTrack>& tracks) const {
27  std::vector<TransientTrack> remaining;
28 
29  return vertices(tracks, remaining, reco::BeamSpot(), false);
30 }
31 
32 std::vector<TransientVertex> ConfigurableTrimmedVertexFinder::vertices(const std::vector<TransientTrack>& tracks,
33  const reco::BeamSpot& spot) const {
34  std::vector<TransientTrack> remaining;
35  return vertices(tracks, remaining, spot, true);
36 }
37 
38 std::vector<TransientVertex> ConfigurableTrimmedVertexFinder::vertices(const std::vector<TransientTrack>& tracks,
39  std::vector<TransientTrack>& unused,
40  const reco::BeamSpot& spot,
41  bool use_spot) const {
44 
45  std::vector<TransientTrack> filtered;
46  for (std::vector<TransientTrack>::const_iterator it = tracks.begin(); it != tracks.end(); it++) {
47  if (theFilter(*it)) {
48  filtered.push_back(*it);
49  } else {
50  unused.push_back(*it);
51  }
52  }
53 
54  std::vector<TransientVertex> all = vertexCandidates(filtered, unused, spot, use_spot);
55 
57 
58  std::vector<TransientVertex> sel = clean(all);
59 
61 
62  return sel;
63 }
64 
65 std::vector<TransientVertex> ConfigurableTrimmedVertexFinder::vertexCandidates(const std::vector<TransientTrack>& tracks,
66  std::vector<TransientTrack>& unused,
67  const reco::BeamSpot& spot,
68  bool use_spot) const {
69  std::vector<TransientVertex> cand;
70 
71  std::vector<TransientTrack> remain = tracks;
72 
73  while (true) {
74  float tkCompCut = (cand.empty() ? theTrackCompatibilityToPV : theTrackCompatibilityToSV);
75 
76  // std::cout << "PVR:compat cut " << tkCompCut << std::endl;
78  // std::cout << "PVCF:compat cut after setting "
79  // << theClusterFinder.trackCompatibilityCut() << std::endl;
80 
81  std::vector<TransientVertex> newVertices;
82  if (cand.empty() && use_spot) {
83  newVertices = theClusterFinder.vertices(remain, spot);
84  } else {
85  newVertices = theClusterFinder.vertices(remain);
86  }
87  if (newVertices.empty())
88  break;
89 
90  analyseClusterFinder(newVertices, remain);
91 
92  for (std::vector<TransientVertex>::const_iterator iv = newVertices.begin(); iv != newVertices.end(); iv++) {
93  if (iv->originalTracks().size() > 1) {
94  cand.push_back(*iv);
95  } else {
96  // candidate has too few tracks - get them back into the vector
97  for (std::vector<TransientTrack>::const_iterator trk = iv->originalTracks().begin();
98  trk != iv->originalTracks().end();
99  ++trk) {
100  unused.push_back(*trk);
101  }
102  }
103  }
104 
105  // when max number of vertices reached, stop
106  if (theMaxNbOfVertices != 0) {
107  if (cand.size() >= (unsigned int)theMaxNbOfVertices)
108  break;
109  }
110  }
111 
112  for (std::vector<TransientTrack>::const_iterator it = remain.begin(); it != remain.end(); it++) {
113  unused.push_back(*it);
114  }
115 
116  return cand;
117 }
118 
119 std::vector<TransientVertex> ConfigurableTrimmedVertexFinder::clean(
120  const std::vector<TransientVertex>& candidates) const {
121  std::vector<TransientVertex> sel;
122  for (std::vector<TransientVertex>::const_iterator i = candidates.begin(); i != candidates.end(); i++) {
123  if (ChiSquaredProbability((*i).totalChiSquared(), (*i).degreesOfFreedom()) > theVtxFitProbCut) {
124  sel.push_back(*i);
125  }
126  }
127 
128  return sel;
129 }
ConfigurableTrimmedVertexFinder::analyseFoundVertices
virtual void analyseFoundVertices(const std::vector< TransientVertex > &vts) const
Definition: ConfigurableTrimmedVertexFinder.h:95
PDWG_EXOHSCP_cff.tracks
tracks
Definition: PDWG_EXOHSCP_cff.py:28
VertexUpdator< 5 >
mps_fire.i
i
Definition: mps_fire.py:355
ConfigurableTrimmedVertexFinder::analyseClusterFinder
virtual void analyseClusterFinder(const std::vector< TransientVertex > &vts, const std::vector< reco::TransientTrack > &remain) const
Definition: ConfigurableTrimmedVertexFinder.h:90
ChiSquaredProbability
float ChiSquaredProbability(double chiSquared, double nrDOF)
Definition: ChiSquaredProbability.cc:13
TrimmedTrackFilter::setPtCut
void setPtCut(double ptCut)
Definition: TrimmedTrackFilter.h:27
ChiSquaredProbability.h
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
ConfigurableTrimmedVertexFinder::theTrackCompatibilityToPV
float theTrackCompatibilityToPV
Definition: ConfigurableTrimmedVertexFinder.h:116
python.cmstools.all
def all(container)
workaround iterator generators for ROOT classes
Definition: cmstools.py:26
VertexFitter< 5 >
ConfigurableTrimmedVertexFinder::ConfigurableTrimmedVertexFinder
ConfigurableTrimmedVertexFinder(const VertexFitter< 5 > *vf, const VertexUpdator< 5 > *vu, const VertexTrackCompatibilityEstimator< 5 > *ve)
Definition: ConfigurableTrimmedVertexFinder.cc:6
alignCSCRings.s
s
Definition: alignCSCRings.py:92
ConfigurableTrimmedVertexFinder::theVtxFitProbCut
float theVtxFitProbCut
Definition: ConfigurableTrimmedVertexFinder.h:115
reco::BeamSpot
Definition: BeamSpot.h:21
ConfigurableTrimmedVertexFinder::vertexCandidates
std::vector< TransientVertex > vertexCandidates(const std::vector< reco::TransientTrack > &tracks, std::vector< reco::TransientTrack > &unused, const reco::BeamSpot &spot, bool use_spot) const
Definition: ConfigurableTrimmedVertexFinder.cc:65
ConfigurableTrimmedVertexFinder::analyseVertexCandidates
virtual void analyseVertexCandidates(const std::vector< TransientVertex > &vts) const
Definition: ConfigurableTrimmedVertexFinder.h:93
edm::ParameterSet
Definition: ParameterSet.h:36
ConfigurableTrimmedVertexFinder::clean
std::vector< TransientVertex > clean(const std::vector< TransientVertex > &candidates) const
Definition: ConfigurableTrimmedVertexFinder.cc:119
cand
Definition: decayParser.h:34
createfilelist.int
int
Definition: createfilelist.py:10
filtered
static const TGPicture * filtered(bool iBackgroundIsBlack)
Definition: FWCollectionSummaryWidget.cc:99
ConfigurableTrimmedVertexFinder::resetEvent
virtual void resetEvent(const std::vector< reco::TransientTrack > &tracks) const
Definition: ConfigurableTrimmedVertexFinder.h:86
TrimmedVertexFinder::vertices
std::vector< TransientVertex > vertices(std::vector< reco::TransientTrack > &remain) const
Definition: TrimmedVertexFinder.cc:26
VertexTrackCompatibilityEstimator< 5 >
HLT_2018_cff.candidates
candidates
Definition: HLT_2018_cff.py:53513
ConfigurableTrimmedVertexFinder::theFilter
TrimmedTrackFilter theFilter
Definition: ConfigurableTrimmedVertexFinder.h:119
ConfigurableTrimmedVertexFinder::analyseInputTracks
virtual void analyseInputTracks(const std::vector< reco::TransientTrack > &tracks) const
Definition: ConfigurableTrimmedVertexFinder.h:88
ConfigurableTrimmedVertexFinder::theClusterFinder
TrimmedVertexFinder theClusterFinder
Definition: ConfigurableTrimmedVertexFinder.h:114
ConfigurableTrimmedVertexFinder.h
EgammaValidation_Wenu_cff.sel
sel
Definition: EgammaValidation_Wenu_cff.py:33
ConfigurableTrimmedVertexFinder::theTrackCompatibilityToSV
float theTrackCompatibilityToSV
Definition: ConfigurableTrimmedVertexFinder.h:117
TrimmedVertexFinder::setTrackCompatibilityCut
void setTrackCompatibilityCut(float cut)
Definition: TrimmedVertexFinder.h:61
ConfigurableTrimmedVertexFinder::theMaxNbOfVertices
int theMaxNbOfVertices
Definition: ConfigurableTrimmedVertexFinder.h:118
ConfigurableTrimmedVertexFinder::setParameters
void setParameters(const edm::ParameterSet &)
Definition: ConfigurableTrimmedVertexFinder.cc:18
ConfigurableTrimmedVertexFinder::vertices
std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &tracks) const override
Definition: ConfigurableTrimmedVertexFinder.cc:26