CMS 3D CMS Logo

vertexFinder.h
Go to the documentation of this file.
1 #ifndef RecoTracker_PixelVertexFinding_plugins_alpaka_vertexFinder_h
2 #define RecoTracker_PixelVertexFinding_plugins_alpaka_vertexFinder_h
3 
4 #include <cstddef>
5 #include <cstdint>
6 
7 #include <alpaka/alpaka.hpp>
8 
17 
19 
20  using namespace cms::alpakatools;
23 
24  class Init {
25  public:
26  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
27  ALPAKA_FN_ACC void operator()(const TAcc &acc, VtxSoAView pdata, WsSoAView pws) const {
28  pdata.nvFinal() = 0; // initialization
30  }
31  };
32 
33  template <typename TrackerTraits>
34  class Producer {
36 
37  public:
39  bool useDensity,
40  bool useDBSCAN,
41  bool useIterative,
42  bool doSplitting,
43  int iminT, // min number of neighbours to be "core"
44  float ieps, // max absolute distance to cluster
45  float ierrmax, // max error to be "seed"
46  float ichi2max // max normalized distance to cluster
47  )
48  : oneKernel_(oneKernel && !(useDBSCAN || useIterative)),
49  useDensity_(useDensity),
50  useDBSCAN_(useDBSCAN),
51  useIterative_(useIterative),
52  doSplitting_(doSplitting),
53  minT(iminT),
54  eps(ieps),
55  errmax(ierrmax),
56  chi2max(ichi2max) {}
57 
58  ~Producer() = default;
59 
60  ZVertexSoACollection makeAsync(Queue &queue, const TkSoAConstView &tracks_view, float ptMin, float ptMax) const;
61 
62  private:
63  const bool oneKernel_; // run everything (cluster,fit,split,sort) in one kernel. Uses only density clusterizer
64  const bool useDensity_; // use density clusterizer
65  const bool useDBSCAN_; // use DBScan clusterizer
66  const bool useIterative_; // use iterative clusterizer
67  const bool doSplitting_; //run vertex splitting
68 
69  int minT; // min number of neighbours to be "core"
70  float eps; // max absolute distance to cluster
71  float errmax; // max error to be "seed"
72  float chi2max; // max normalized distance to cluster
73  };
74 
75 } // namespace ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder
76 
77 #endif // RecoTracker_PixelVertexFinding_plugins_alpaka_vertexFinder_h
ALPAKA_FN_ACC ALPAKA_FN_INLINE void VtxSoAView & pdata
ALPAKA_FN_ACC ALPAKA_FN_INLINE void VtxSoAView WsSoAView int float eps
std::conditional_t< std::is_same_v< Device, alpaka::DevCpu >, ZVertexHost, ZVertexDevice< Device > > ZVertexSoACollection
uint32_t const *__restrict__ TkSoAView< TrackerTraits > tracks_view
Producer(bool oneKernel, bool useDensity, bool useDBSCAN, bool useIterative, bool doSplitting, int iminT, float ieps, float ierrmax, float ichi2max)
Definition: vertexFinder.h:38
ALPAKA_FN_ACC ALPAKA_FN_INLINE void VtxSoAView WsSoAView int float float errmax
int init
Definition: HydjetWrapper.h:66
constexpr float ptMin
ALPAKA_FN_ACC void operator()(const TAcc &acc, VtxSoAView pdata, WsSoAView pws) const
Definition: vertexFinder.h:27
ZVertexSoAHeterogeneousLayout<>::View ZVertexSoAView
reco::TrackSoAConstView< TrackerTraits > TkSoAConstView
Definition: vertexFinder.h:35
ALPAKA_FN_ACC ALPAKA_FN_INLINE void VtxSoAView WsSoAView & pws
ALPAKA_FN_ACC ALPAKA_FN_INLINE void VtxSoAView WsSoAView int minT
PixelVertexWSSoALayout<>::View PixelVertexWorkSpaceSoAView
::vertexFinder::PixelVertexWorkSpaceSoAView WsSoAView
typename reco::TrackSoA< TrackerTraits >::template Layout<>::ConstView TrackSoAConstView
Definition: TracksSoA.h:47
ALPAKA_FN_ACC ALPAKA_FN_INLINE void VtxSoAView WsSoAView int float float float chi2max