CMS 3D CMS Logo

HelixFitOnGPU.h
Go to the documentation of this file.
1 #ifndef RecoTracker_PixelSeeding_plugins_HelixFitOnGPU_h
2 #define RecoTracker_PixelSeeding_plugins_HelixFitOnGPU_h
3 
8 
9 #include "CAStructures.h"
10 
11 namespace riemannFit {
12  // in case of memory issue can be made smaller
13  constexpr uint32_t maxNumberOfConcurrentFits = 32 * 1024;
15  using Matrix3x4d = Eigen::Matrix<double, 3, 4>;
16  using Map3x4d = Eigen::Map<Matrix3x4d, 0, Eigen::Stride<3 * stride, stride> >;
17  using Matrix6x4f = Eigen::Matrix<float, 6, 4>;
18  using Map6x4f = Eigen::Map<Matrix6x4f, 0, Eigen::Stride<6 * stride, stride> >;
19 
20  // hits
21  template <int N>
22  using Matrix3xNd = Eigen::Matrix<double, 3, N>;
23  template <int N>
24  using Map3xNd = Eigen::Map<Matrix3xNd<N>, 0, Eigen::Stride<3 * stride, stride> >;
25  // errors
26  template <int N>
27  using Matrix6xNf = Eigen::Matrix<float, 6, N>;
28  template <int N>
29  using Map6xNf = Eigen::Map<Matrix6xNf<N>, 0, Eigen::Stride<6 * stride, stride> >;
30  // fast fit
31  using Map4d = Eigen::Map<Vector4d, 0, Eigen::InnerStride<stride> >;
32 
33  template <auto Start, auto End, auto Inc, class F> //a compile-time bounded for loop
34  constexpr void rolling_fits(F &&f) {
35  if constexpr (Start < End) {
36  f(std::integral_constant<decltype(Start), Start>());
37  rolling_fits<Start + Inc, End, Inc>(f);
38  }
39  }
40 
41 } // namespace riemannFit
42 
43 template <typename TrackerTraits>
45 public:
47 
50 
53 
55 
56  explicit HelixFitOnGPU(float bf, bool fitNas4) : bField_(bf), fitNas4_(fitNas4) {}
58 
59  void setBField(double bField) { bField_ = bField; }
60  void launchRiemannKernels(const HitConstView &hv, uint32_t nhits, uint32_t maxNumberOfTuples, cudaStream_t cudaStream);
62  uint32_t nhits,
63  uint32_t maxNumberOfTuples,
64  cudaStream_t cudaStream);
65 
66  void launchRiemannKernelsOnCPU(const HitConstView &hv, uint32_t nhits, uint32_t maxNumberOfTuples);
67  void launchBrokenLineKernelsOnCPU(const HitConstView &hv, uint32_t nhits, uint32_t maxNumberOfTuples);
68 
69  void allocateOnGPU(TupleMultiplicity const *tupleMultiplicity, OutputSoAView &helix_fit_results);
70  void deallocateOnGPU();
71 
72 private:
74 
75  // fowarded
76  Tuples const *tuples_ = nullptr;
79  float bField_;
80 
81  const bool fitNas4_;
82 };
83 
84 #endif // RecoTracker_PixelSeeding_plugins_HelixFitOnGPU_h
Eigen::Matrix< double, 3, 4 > Matrix3x4d
Definition: HelixFit.h:23
Eigen::Map< Matrix3x4d, 0, Eigen::Stride< 3 *stride, stride > > Map3x4d
Definition: HelixFit.h:24
typename TrackSoA< TrackerTraits >::HitContainer Tuples
Definition: HelixFitOnGPU.h:51
Eigen::Matrix< float, 6, 4 > Matrix6x4f
Definition: HelixFit.h:25
TupleMultiplicity< TrackerTraits > const *__restrict__ tupleMultiplicity
constexpr uint32_t stride
Definition: HelixFit.h:22
const bool fitNas4_
Definition: HelixFitOnGPU.h:81
Eigen::Matrix< float, 6, N > Matrix6xNf
Definition: HelixFit.h:35
constexpr void rolling_fits(F &&f)
Definition: HelixFit.h:42
Eigen::Map< Vector4d, 0, Eigen::InnerStride< stride > > Map4d
Definition: HelixFit.h:39
Eigen::Map< Matrix3xNd< N >, 0, Eigen::Stride< 3 *stride, stride > > Map3xNd
Definition: HelixFit.h:32
void launchRiemannKernelsOnCPU(const HitConstView &hv, uint32_t nhits, uint32_t maxNumberOfTuples)
TrackingRecHitSoAView< TrackerTraits > HitView
Definition: HelixFitOnGPU.h:48
Eigen::Map< Matrix6x4f, 0, Eigen::Stride< 6 *stride, stride > > Map6x4f
Definition: HelixFit.h:26
static constexpr uint32_t maxNumberOfConcurrentFits_
Definition: HelixFitOnGPU.h:73
void setBField(double bField)
Definition: HelixFitOnGPU.h:59
HelixFitOnGPU(float bf, bool fitNas4)
Definition: HelixFitOnGPU.h:56
double f[11][100]
void launchRiemannKernels(const HitConstView &hv, uint32_t nhits, uint32_t maxNumberOfTuples, cudaStream_t cudaStream)
Eigen::Matrix< double, 3, N > Matrix3xNd
Definition: HelixFit.h:30
Eigen::Map< Matrix6xNf< N >, 0, Eigen::Stride< 6 *stride, stride > > Map6xNf
Definition: HelixFit.h:37
void launchBrokenLineKernelsOnCPU(const HitConstView &hv, uint32_t nhits, uint32_t maxNumberOfTuples)
TrackSoAView< TrackerTraits > OutputSoAView
Definition: HelixFitOnGPU.h:52
TupleMultiplicity const * tupleMultiplicity_
Definition: HelixFitOnGPU.h:77
typename TrackingRecHitSoA< TrackerTraits >::template TrackingRecHitSoALayout<>::ConstView TrackingRecHitSoAConstView
void allocateOnGPU(TupleMultiplicity const *tupleMultiplicity, OutputSoAView &helix_fit_results)
Definition: HelixFitOnGPU.cc:5
Tuples const * tuples_
Definition: HelixFitOnGPU.h:76
typename TrackSoA< TrackerTraits >::template TrackSoALayout<>::View TrackSoAView
TrackingRecHitSoAConstView< TrackerTraits > HitConstView
Definition: HelixFitOnGPU.h:49
void deallocateOnGPU()
typename TrackingRecHitSoA< TrackerTraits >::template TrackingRecHitSoALayout<>::View TrackingRecHitSoAView
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
Definition: blowfish.cc:163
OutputSoAView outputSoa_
Definition: HelixFitOnGPU.h:78
void launchBrokenLineKernels(const HitConstView &hv, uint32_t nhits, uint32_t maxNumberOfTuples, cudaStream_t cudaStream)
reco::TrackSoAView< TrackerTraits > OutputSoAView
constexpr uint32_t maxNumberOfConcurrentFits
Definition: HelixFit.h:21