HeterogeneousCore
CUDAUtilities
interface
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
10
#include "
HeterogeneousCore/CUDAUtilities/interface/cudaCompat.h
"
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
21
__host__
__device__
constexpr
Scalar
&
operator()
(int32_t
i
) {
return
data_
[
i
]; }
22
__device__
constexpr
const
Scalar
operator()
(int32_t
i
)
const
{
return
__ldg
(
data_
+
i
); }
23
__host__
__device__
constexpr
Scalar
&
operator[]
(int32_t
i
) {
return
data_
[
i
]; }
24
__device__
constexpr
const
Scalar
operator[]
(int32_t
i
)
const
{
return
__ldg
(
data_
+
i
); }
25
26
__host__
__device__
constexpr
Scalar
*
data
() {
return
data_
; }
27
__host__
__device__
constexpr
Scalar
const
*
data
()
const
{
return
data_
; }
28
29
private
:
30
Scalar
data_
[
S
];
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
42
__host__
__device__
constexpr
Map
operator()
(int32_t
i
) {
return
Map
(
data_
+
i
); }
43
__host__
__device__
constexpr
CMap
operator()
(int32_t
i
)
const
{
return
CMap
(
data_
+
i
); }
44
__host__
__device__
constexpr
Map
operator[]
(int32_t
i
) {
return
Map
(
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
eigenSoA::MatrixSoA::Map
Eigen::Map< M, 0, Eigen::Stride< M::RowsAtCompileTime *S, S > > Map
Definition:
eigenSoA.h:39
mps_fire.i
i
Definition:
mps_fire.py:428
eigenSoA::isPowerOf2
constexpr bool isPowerOf2(int32_t v)
Definition:
eigenSoA.h:14
eigenSoA::ScalarSoA::data
constexpr Scalar * data()
Definition:
eigenSoA.h:26
eigenSoA::ScalarSoA
Definition:
eigenSoA.h:17
__device__
#define __device__
Definition:
cudaCompat.h:92
findQualityFiles.v
v
Definition:
findQualityFiles.py:179
align::Scalar
double Scalar
Definition:
Definitions.h:25
eigenSoA::MatrixSoA::CMap
Eigen::Map< const M, 0, Eigen::Stride< M::RowsAtCompileTime *S, S > > CMap
Definition:
eigenSoA.h:40
eigenSoA::MatrixSoA::operator()
constexpr CMap operator()(int32_t i) const
Definition:
eigenSoA.h:43
eigenSoA::ScalarSoA::operator()
constexpr Scalar & operator()(int32_t i)
Definition:
eigenSoA.h:21
eigenSoA::ScalarSoA::data_
Scalar data_[S]
Definition:
eigenSoA.h:30
eigenSoA::MatrixSoA::operator()
constexpr Map operator()(int32_t i)
Definition:
eigenSoA.h:42
eigenSoA::ScalarSoA::data
constexpr Scalar const * data() const
Definition:
eigenSoA.h:27
S
double S(const TLorentzVector &, const TLorentzVector &)
Definition:
Particle.cc:97
eigenSoA::MatrixSoA::operator[]
constexpr Map operator[](int32_t i)
Definition:
eigenSoA.h:44
__host__
#define __host__
Definition:
cudaCompat.h:91
eigenSoA::ScalarSoA::operator[]
constexpr const Scalar operator[](int32_t i) const
Definition:
eigenSoA.h:24
eigenSoA::ScalarSoA::operator[]
constexpr Scalar & operator[](int32_t i)
Definition:
eigenSoA.h:23
eigenSoA::ScalarSoA::operator()
constexpr const Scalar operator()(int32_t i) const
Definition:
eigenSoA.h:22
eigenSoA::MatrixSoA
Definition:
eigenSoA.h:36
eigenSoA::MatrixSoA::operator[]
constexpr CMap operator[](int32_t i) const
Definition:
eigenSoA.h:45
eigenSoA::MatrixSoA::data_
Scalar data_[S *M::RowsAtCompileTime *M::ColsAtCompileTime]
Definition:
eigenSoA.h:48
T
long double T
Definition:
Basic3DVectorLD.h:48
eigenSoA::ScalarSoA::Scalar
T Scalar
Definition:
eigenSoA.h:19
S
Definition:
CSCDBL1TPParametersExtended.h:16
cudaCompat.h
eigenSoA
Definition:
eigenSoA.h:12
cms::cudacompat::__ldg
T __ldg(T const *x)
Definition:
cudaCompat.h:77
eigenSoA::MatrixSoA::Scalar
typename M::Scalar Scalar
Definition:
eigenSoA.h:38
Generated for CMSSW Reference Manual by
1.8.16