CMS 3D CMS Logo

eigenSoA.h
Go to the documentation of this file.
1 #ifndef HeterogeneousCore_CUDAUtilities_interface_eigenSoA_h
2 #define HeterogeneousCore_CUDAUtilities_interface_eigenSoA_h
3 
4 #include <algorithm>
5 #include <cmath>
6 #include <cstdint>
7 
8 #include <Eigen/Core>
9 
11 
12 namespace eigenSoA {
13 
14  constexpr bool isPowerOf2(int32_t v) { return v && !(v & (v - 1)); }
15 
16  template <typename T, int S>
17  class alignas(128) ScalarSoA {
18  public:
19  using Scalar = T;
20 
22  __device__ constexpr const Scalar operator()(int32_t i) const { return __ldg(data_ + i); }
24  __device__ constexpr const Scalar operator[](int32_t i) const { return __ldg(data_ + i); }
25 
27  __host__ __device__ constexpr Scalar const* data() const { return data_; }
28 
29  private:
31  static_assert(isPowerOf2(S), "SoA stride not a power of 2");
32  static_assert(sizeof(data_) % 128 == 0, "SoA size not a multiple of 128");
33  };
34 
35  template <typename M, int S>
36  class alignas(128) MatrixSoA {
37  public:
38  using Scalar = typename M::Scalar;
39  using Map = Eigen::Map<M, 0, Eigen::Stride<M::RowsAtCompileTime * S, S> >;
40  using CMap = Eigen::Map<const M, 0, Eigen::Stride<M::RowsAtCompileTime * S, S> >;
41 
43  __host__ __device__ constexpr CMap operator()(int32_t i) const { return CMap(data_ + i); }
45  __host__ __device__ constexpr CMap operator[](int32_t i) const { return CMap(data_ + i); }
46 
47  private:
48  Scalar data_[S * M::RowsAtCompileTime * M::ColsAtCompileTime];
49  static_assert(isPowerOf2(S), "SoA stride not a power of 2");
50  static_assert(sizeof(data_) % 128 == 0, "SoA size not a multiple of 128");
51  };
52 
53 } // namespace eigenSoA
54 
55 #endif // HeterogeneousCore_CUDAUtilities_interface_eigenSoA_h
double Scalar
Definition: Definitions.h:25
#define __host__
constexpr bool isPowerOf2(int32_t v)
Definition: eigenSoA.h:14
Scalar data_[S]
Definition: eigenSoA.h:30
Eigen::Map< const M, 0, Eigen::Stride< M::RowsAtCompileTime *S, S > > CMap
Definition: eigenSoA.h:40
__device__ constexpr const Scalar operator()(int32_t i) const
Definition: eigenSoA.h:22
__host__ __device__ constexpr Map operator[](int32_t i)
Definition: eigenSoA.h:44
__host__ __device__ constexpr Map operator()(int32_t i)
Definition: eigenSoA.h:42
__host__ __device__ constexpr Scalar & operator[](int32_t i)
Definition: eigenSoA.h:23
__device__ constexpr const Scalar operator[](int32_t i) const
Definition: eigenSoA.h:24
Scalar data_[S *M::RowsAtCompileTime *M::ColsAtCompileTime]
Definition: eigenSoA.h:48
__host__ __device__ constexpr Scalar & operator()(int32_t i)
Definition: eigenSoA.h:21
typename M::Scalar Scalar
Definition: eigenSoA.h:38
__host__ __device__ constexpr CMap operator[](int32_t i) const
Definition: eigenSoA.h:45
T __ldg(T const *x)
Definition: cudaCompat.h:137
__host__ __device__ constexpr Scalar * data()
Definition: eigenSoA.h:26
__host__ __device__ constexpr CMap operator()(int32_t i) const
Definition: eigenSoA.h:43
__host__ __device__ constexpr Scalar const * data() const
Definition: eigenSoA.h:27
long double T
#define __device__
Eigen::Map< M, 0, Eigen::Stride< M::RowsAtCompileTime *S, S > > Map
Definition: eigenSoA.h:39