CMS 3D CMS Logo

VertexFinder.h
Go to the documentation of this file.
1 #ifndef __L1Trigger_VertexFinder_VertexFinder_h__
2 #define __L1Trigger_VertexFinder_VertexFinder_h__
3 
11 
12 #include <algorithm>
13 #include <iterator>
14 #include <vector>
15 
16 namespace l1tVertexFinder {
17 
18  // Returns the number of bits needed to represent and integer decimal
19  static constexpr unsigned BitsToRepresent(unsigned x) { return x < 2 ? 1 : 1 + BitsToRepresent(x >> 1); }
20 
21  typedef std::vector<L1Track> FitTrackCollection;
22  typedef std::vector<RecoVertex<>> RecoVertexCollection;
23 
24  class VertexFinder {
25  public:
29  settings_ = &settings;
30  }
32 
34  struct SortTracksByZ0 {
35  inline bool operator()(const L1Track track0, const L1Track track1) { return (track0.z0() < track1.z0()); }
36  };
37 
38  struct SortTracksByPt {
39  inline bool operator()(const L1Track track0, const L1Track track1) {
40  return (std::abs(track0.pt()) > std::abs(track1.pt()));
41  }
42  };
43 
45 
47  const FitTrackCollection& fitTracks() const { return fitTracks_; }
49  unsigned int iterationsPerTrack() const { return double(iterations_) / double(fitTracks_.size()); }
51  unsigned int numInputTracks() const { return fitTracks_.size(); }
53  unsigned int numIterations() const { return iterations_; }
55  unsigned int numVertices() const { return vertices_.size(); }
57  unsigned int numVerticesEmulation() const { return verticesEmulation_.size(); }
59  template <typename T>
60  T PrimaryVertex() const {
62  return vertices_[pv_index_];
63  else if ((settings_->vx_precision() == Precision::Emulation) && (pv_index_ < vertices_.size()))
65  else {
66  edm::LogWarning("VertexFinder") << "PrimaryVertex::No primary vertex has been found.";
67  return RecoVertex<>();
68  }
69  }
71  unsigned int primaryVertexId() const { return pv_index_; }
73  const std::vector<RecoVertex<>>& vertices() const { return vertices_; }
76 
78  void findPrimaryVertex();
80  void associatePrimaryVertex(double trueZ0);
82  void GapClustering();
84  float maxDistance(RecoVertex<> cluster0, RecoVertex<> cluster1);
86  float minDistance(RecoVertex<> cluster0, RecoVertex<> cluster1);
88  float meanDistance(RecoVertex<> cluster0, RecoVertex<> cluster1);
90  float centralDistance(RecoVertex<> cluster0, RecoVertex<> cluster1);
98  void fastHisto(const TrackerTopology* tTopo);
100  void fastHistoEmulation();
101 
103  void sortVerticesInPt();
105  void sortVerticesInZ0();
106 
108  template <class data_type, typename stream_type = std::ostream>
109  void printHistogram(stream_type& stream,
110  std::vector<data_type> data,
111  int width = 80,
112  int minimum = 0,
113  int maximum = -1,
114  std::string title = "",
115  std::string color = "");
116 
117  template <typename ForwardIterator, typename T>
118  void strided_iota(ForwardIterator first, ForwardIterator last, T value, T stride) {
119  while (first != last) {
120  *first++ = value;
121  value += stride;
122  }
123  }
124 
126 
129  const std::vector<float>& bin_centers,
130  const std::vector<unsigned int>& counts);
132  void DBSCAN();
134  void HPV();
136  void Kmeans();
138  void PVR();
139 
140  private:
144  unsigned int numMatchedVertices_;
146  unsigned int pv_index_;
147  unsigned int iterations_;
148  };
149 
150 } // end namespace l1tVertexFinder
151 
152 #endif
bool operator()(const L1Track track0, const L1Track track1)
Definition: VertexFinder.h:35
const FitTrackCollection & fitTracks() const
Accessors.
Definition: VertexFinder.h:47
FitTrackCollection fitTracks_
Definition: VertexFinder.h:145
void fastHistoLooseAssociation()
High pT Vertex Algorithm.
float maxDistance(RecoVertex<> cluster0, RecoVertex<> cluster1)
Find maximum distance in two clusters of tracks.
Definition: VertexFinder.cc:83
void GapClustering()
Gap Clustering Algorithm.
Definition: VertexFinder.cc:65
float meanDistance(RecoVertex<> cluster0, RecoVertex<> cluster1)
Find average distance in two clusters of tracks.
bool operator()(const L1Track track0, const L1Track track1)
Definition: VertexFinder.h:39
float pt() const
Definition: L1Track.h:24
std::vector< L1Track > FitTrackCollection
Definition: VertexFinder.h:21
float centralDistance(RecoVertex<> cluster0, RecoVertex<> cluster1)
Find distance between centres of two clusters.
unsigned int numVerticesEmulation() const
Number of emulation vertices.
Definition: VertexFinder.h:57
uint32_t T const *__restrict__ uint32_t const *__restrict__ int32_t int Histo::index_type cudaStream_t stream
void agglomerativeHierarchicalClustering()
Simple Merge Algorithm.
T PrimaryVertex() const
Reconstructed primary vertex.
Definition: VertexFinder.h:60
void adaptiveVertexReconstruction()
Adaptive Vertex Reconstruction algorithm.
void findPrimaryVertex()
Find the primary vertex.
void fastHisto(const TrackerTopology *tTopo)
Histogramming algorithm.
static constexpr unsigned BitsToRepresent(unsigned x)
Definition: VertexFinder.h:19
const std::vector< RecoVertex<> > & vertices() const
Returns the z positions of the reconstructed primary vertices.
Definition: VertexFinder.h:73
Simple wrapper class for TTTrack.
Definition: L1Track.h:17
void printHistogram(stream_type &stream, std::vector< data_type > data, int width=80, int minimum=0, int maximum=-1, std::string title="", std::string color="")
Print an ASCII histogram.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< VertexWord > VertexWordCollection
Definition: VertexWord.h:195
Definition: value.py:1
VertexFinder(FitTrackCollection &fitTracks, const AlgoSettings &settings)
Constructor and destructor.
Definition: VertexFinder.h:27
RecoVertexCollection vertices_
Definition: VertexFinder.h:142
unsigned int numInputTracks() const
Storage for tracks out of the L1 Track finder.
Definition: VertexFinder.h:51
float minDistance(RecoVertex<> cluster0, RecoVertex<> cluster1)
Find minimum distance in two clusters of tracks.
Definition: VertexFinder.cc:96
void sortVerticesInPt()
Sort vertices in pT.
void PVR()
Find maximum distance in two clusters of tracks.
void strided_iota(ForwardIterator first, ForwardIterator last, T value, T stride)
Definition: VertexFinder.h:118
const AlgoSettings * settings_
Definition: VertexFinder.h:141
Precision vx_precision() const
Definition: AlgoSettings.h:35
void fastHistoEmulation()
Histogramming algorithm (emulation)
float z0() const
Definition: L1Track.h:25
void HPV()
High pT Vertex Algorithm.
std::vector< RecoVertex<> > RecoVertexCollection
Definition: VertexFinder.h:22
void sortVerticesInZ0()
Sort vertices in z.
l1t::VertexWordCollection verticesEmulation_
Definition: VertexFinder.h:143
unsigned int numIterations() const
Number of iterations.
Definition: VertexFinder.h:53
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:79
void Kmeans()
Kmeans Algorithm.
void associatePrimaryVertex(double trueZ0)
Associate the primary vertex with the real one.
const l1t::VertexWordCollection & verticesEmulation() const
Returns the emulation primary vertices.
Definition: VertexFinder.h:75
unsigned int primaryVertexId() const
Reconstructed Primary Vertex Id.
Definition: VertexFinder.h:71
Log< level::Warning, false > LogWarning
unsigned int numVertices() const
Number of reconstructed vertices.
Definition: VertexFinder.h:55
long double T
unsigned int iterationsPerTrack() const
Number of iterations.
Definition: VertexFinder.h:49
void DBSCAN()
DBSCAN algorithm.
void computeAndSetVertexParameters(RecoVertex<> &vertex, const std::vector< float > &bin_centers, const std::vector< unsigned int > &counts)
Vertexing algorithms.
Definition: VertexFinder.cc:7