CMS 3D CMS Logo

vertexFinder.dev.cc
Go to the documentation of this file.
1 #include <alpaka/alpaka.hpp>
2 
8 
9 #include "vertexFinder.h"
10 #include "clusterTracksDBSCAN.h"
11 #include "clusterTracksIterative.h"
12 #include "clusterTracksByDensity.h"
13 #include "fitVertices.h"
14 #include "sortByPt2.h"
15 #include "splitVertices.h"
16 
17 #undef PIXVERTEX_DEBUG_PRODUCE
19  namespace vertexFinder {
20  using namespace cms::alpakatools;
21  // reject outlier tracks that contribute more than this to the chi2 of the vertex fit
23  constexpr float maxChi2ForFinalFit = 5000.f;
24 
25  // split vertices with a chi2/NDoF greater than this
27 
28  template <typename TrackerTraits>
29  class LoadTracks {
30  public:
31  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
32  ALPAKA_FN_ACC void operator()(const TAcc& acc,
35  WsSoAView pws,
36  float ptMin,
37  float ptMax) const {
38  auto const* quality = tracks_view.quality();
40 
41  for (auto idx : cms::alpakatools::uniform_elements(acc, tracks_view.nTracks())) {
42  [[maybe_unused]] auto nHits = helper::nHits(tracks_view, idx);
44 
45  // initialize soa...
46  soa[idx].idv() = -1;
47 
49  continue; // no triplets
51  continue;
52 
53  auto pt = tracks_view[idx].pt();
54 
55  if (pt < ptMin)
56  continue;
57 
58  // clamp pt
59  pt = std::min<float>(pt, ptMax);
60 
61  auto& data = pws;
62  auto it = alpaka::atomicAdd(acc, &data.ntrks(), 1u, alpaka::hierarchy::Blocks{});
63  data[it].itrk() = idx;
64  data[it].zt() = reco::zip(tracks_view, idx);
65  data[it].ezt2() = tracks_view[idx].covariance()(14);
66  data[it].ptt2() = pt * pt;
67  }
68  }
69  };
70 // #define THREE_KERNELS
71 #ifndef THREE_KERNELS
73  public:
74  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
75  ALPAKA_FN_ACC void operator()(const TAcc& acc,
77  WsSoAView pws,
78  bool doSplit,
79  int minT, // min number of neighbours to be "seed"
80  float eps, // max absolute distance to cluster
81  float errmax, // max error to be "seed"
82  float chi2max // max normalized distance to cluster,
83  ) const {
84  clusterTracksByDensity(acc, pdata, pws, minT, eps, errmax, chi2max);
85  alpaka::syncBlockThreads(acc);
87  alpaka::syncBlockThreads(acc);
88  if (doSplit) {
90  alpaka::syncBlockThreads(acc);
92  alpaka::syncBlockThreads(acc);
93  }
94  sortByPt2(acc, pdata, pws);
95  }
96  };
97 #else
98  class VertexFinderKernel1 {
99  public:
100  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
101  ALPAKA_FN_ACC void operator()(const TAcc& acc,
103  WsSoAView pws,
104  int minT, // min number of neighbours to be "seed"
105  float eps, // max absolute distance to cluster
106  float errmax, // max error to be "seed"
107  float chi2max // max normalized distance to cluster,
108  ) const {
109  clusterTracksByDensity(pdata, pws, minT, eps, errmax, chi2max);
110  alpaka::syncBlockThreads(acc);
112  }
113  };
114  class VertexFinderKernel2 {
115  public:
116  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
117  ALPAKA_FN_ACC void operator()(const TAcc& acc, VtxSoAView pdata, WsSoAView pws) const {
119  alpaka::syncBlockThreads(acc);
120  sortByPt2(pdata, pws);
121  }
122  };
123 #endif
124 
125  template <typename TrackerTraits>
128  float ptMin,
129  float ptMax) const {
130 #ifdef PIXVERTEX_DEBUG_PRODUCE
131  std::cout << "producing Vertices on GPU" << std::endl;
132 #endif // PIXVERTEX_DEBUG_PRODUCE
134 
135  auto soa = vertices.view();
136 
138 
139  // Initialize
140  const auto initWorkDiv = cms::alpakatools::make_workdiv<Acc1D>(1, 1);
141  alpaka::exec<Acc1D>(queue, initWorkDiv, Init{}, soa, ws_d.view());
142 
143  // Load Tracks
144  const uint32_t blockSize = 128;
145  const uint32_t numberOfBlocks =
146  cms::alpakatools::divide_up_by(tracks_view.metadata().size() + blockSize - 1, blockSize);
147  const auto loadTracksWorkDiv = cms::alpakatools::make_workdiv<Acc1D>(numberOfBlocks, blockSize);
148  alpaka::exec<Acc1D>(
149  queue, loadTracksWorkDiv, LoadTracks<TrackerTraits>{}, tracks_view, soa, ws_d.view(), ptMin, ptMax);
150 
151  // Running too many thread lead to problems when printf is enabled.
152  const auto finderSorterWorkDiv = cms::alpakatools::make_workdiv<Acc1D>(1, 1024 - 128);
153  const auto splitterFitterWorkDiv = cms::alpakatools::make_workdiv<Acc1D>(1024, 128);
154 
155  if (oneKernel_) {
156  // implemented only for density clustesrs
157 #ifndef THREE_KERNELS
158  alpaka::exec<Acc1D>(queue,
159  finderSorterWorkDiv,
160  VertexFinderOneKernel{},
161  soa,
162  ws_d.view(),
163  doSplitting_,
164  minT,
165  eps,
166  errmax,
167  chi2max);
168 #else
169  alpaka::exec<Acc1D>(
170  queue, finderSorterWorkDiv, VertexFinderOneKernel{}, soa, ws_d.view(), minT, eps, errmax, chi2max);
171 
172  // one block per vertex...
173  if (doSplitting_)
174  alpaka::exec<Acc1D>(queue, splitterFitterWorkDiv, SplitVerticesKernel{}, soa, ws_d.view(), maxChi2ForSplit);
175  alpaka::exec<Acc1D>(queue, finderSorterWorkDiv{}, soa, ws_d.view());
176 #endif
177  } else { // five kernels
178  if (useDensity_) {
179  alpaka::exec<Acc1D>(
180  queue, finderSorterWorkDiv, ClusterTracksByDensityKernel{}, soa, ws_d.view(), minT, eps, errmax, chi2max);
181 
182  } else if (useDBSCAN_) {
183  alpaka::exec<Acc1D>(
184  queue, finderSorterWorkDiv, ClusterTracksDBSCAN{}, soa, ws_d.view(), minT, eps, errmax, chi2max);
185  } else if (useIterative_) {
186  alpaka::exec<Acc1D>(
187  queue, finderSorterWorkDiv, ClusterTracksIterative{}, soa, ws_d.view(), minT, eps, errmax, chi2max);
188  }
189  alpaka::exec<Acc1D>(queue, finderSorterWorkDiv, FitVerticesKernel{}, soa, ws_d.view(), maxChi2ForFirstFit);
190 
191  // one block per vertex...
192  if (doSplitting_) {
193  alpaka::exec<Acc1D>(queue, splitterFitterWorkDiv, SplitVerticesKernel{}, soa, ws_d.view(), maxChi2ForSplit);
194 
195  alpaka::exec<Acc1D>(queue, finderSorterWorkDiv, FitVerticesKernel{}, soa, ws_d.view(), maxChi2ForFinalFit);
196  }
197  alpaka::exec<Acc1D>(queue, finderSorterWorkDiv, SortByPt2Kernel{}, soa, ws_d.view());
198  }
199 
200  return vertices;
201  }
202 
203  template class Producer<pixelTopology::Phase1>;
204  template class Producer<pixelTopology::Phase2>;
205  template class Producer<pixelTopology::HIonPhase1>;
206  } // namespace vertexFinder
207 } // namespace ALPAKA_ACCELERATOR_NAMESPACE
ALPAKA_FN_ACC auto uniform_elements(TAcc const &acc, TArgs... args)
Definition: workdivision.h:311
Definition: helper.py:1
ALPAKA_FN_ACC ALPAKA_FN_INLINE void VtxSoAView & pdata
ALPAKA_FN_ACC ALPAKA_FN_INLINE void VtxSoAView WsSoAView int float eps
constexpr Idx divide_up_by(Idx value, Idx divisor)
Definition: workdivision.h:20
ALPAKA_FN_ACC void operator()(const TAcc &acc, reco::TrackSoAConstView< TrackerTraits > tracks_view, VtxSoAView soa, WsSoAView pws, float ptMin, float ptMax) const
std::conditional_t< std::is_same_v< Device, alpaka::DevCpu >, ZVertexHost, ZVertexDevice< Device > > ZVertexSoACollection
uint32_t const *__restrict__ TkSoAView< TrackerTraits > tracks_view
ALPAKA_FN_ACC ALPAKA_FN_INLINE void VtxSoAView WsSoAView int float float errmax
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float zip(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:90
constexpr uint32_t MAXTRACKS
fitVertices(pdata, pws, maxChi2ForFirstFit)
constexpr float ptMin
string quality
ALPAKA_FN_ACC ALPAKA_FN_INLINE void VtxSoAView WsSoAView & pws
std::vector< Block > Blocks
Definition: Block.h:99
ALPAKA_FN_ACC ALPAKA_FN_INLINE void VtxSoAView WsSoAView int minT
splitVertices(pdata, pws, maxChi2ForSplit)
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr bool isTriplet(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:95
ALPAKA_FN_ACC void operator()(const TAcc &acc, VtxSoAView pdata, WsSoAView pws, bool doSplit, int minT, float eps, float errmax, float chi2max) const
ALPAKA_FN_ACC ALPAKA_FN_INLINE void sortByPt2(const TAcc &acc, VtxSoAView &data, WsSoAView &ws)
Definition: sortByPt2.h:26
::vertexFinder::PixelVertexWorkSpaceSoAView WsSoAView
typename reco::TrackSoA< TrackerTraits >::template Layout<>::ConstView TrackSoAConstView
Definition: TracksSoA.h:47
PortableCollection<::vertexFinder::PixelVertexWSSoALayout<> > PixelVertexWorkSpaceSoADevice
ALPAKA_FN_ACC ALPAKA_FN_INLINE void VtxSoAView WsSoAView int float float float chi2max
ZVertexSoACollection makeAsync(Queue &queue, const TkSoAConstView &tracks_view, float ptMin, float ptMax) const
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits
T1 atomicAdd(T1 *a, T2 b)
Definition: cudaCompat.h:61