CMS 3D CMS Logo

TrajectoryStateSoAT.h
Go to the documentation of this file.
1 #ifndef CUDADataFormats_Track_TrajectoryStateSOAT_H
2 #define CUDADataFormats_Track_TrajectoryStateSOAT_H
3 
4 #include <Eigen/Dense>
6 
7 template <int32_t S>
9  using Vector5f = Eigen::Matrix<float, 5, 1>;
10  using Vector15f = Eigen::Matrix<float, 15, 1>;
11 
12  using Vector5d = Eigen::Matrix<double, 5, 1>;
13  using Matrix5d = Eigen::Matrix<double, 5, 5>;
14 
15  static constexpr int32_t stride() { return S; }
16 
19 
20  template <typename V3, typename M3, typename V2, typename M2>
22  V3 const& cp, M3 const& ccov, V2 const& lp, M2 const& lcov, float b, int32_t i) {
23  state(i) << cp.template cast<float>(), lp.template cast<float>();
24  state(i)(2) *= b;
25  auto cov = covariance(i);
26  cov(0) = ccov(0, 0);
27  cov(1) = ccov(0, 1);
28  cov(2) = b * float(ccov(0, 2));
29  cov(4) = cov(3) = 0;
30  cov(5) = ccov(1, 1);
31  cov(6) = b * float(ccov(1, 2));
32  cov(8) = cov(7) = 0;
33  cov(9) = b * b * float(ccov(2, 2));
34  cov(11) = cov(10) = 0;
35  cov(12) = lcov(0, 0);
36  cov(13) = lcov(0, 1);
37  cov(14) = lcov(1, 1);
38  }
39 
40  template <typename V5, typename M5>
41  __host__ __device__ inline void copyFromDense(V5 const& v, M5 const& cov, int32_t i) {
42  state(i) = v.template cast<float>();
43  for (int j = 0, ind = 0; j < 5; ++j)
44  for (auto k = j; k < 5; ++k)
45  covariance(i)(ind++) = cov(j, k);
46  }
47 
48  template <typename V5, typename M5>
49  __host__ __device__ inline void copyToDense(V5& v, M5& cov, int32_t i) const {
50  v = state(i).template cast<typename V5::Scalar>();
51  for (int j = 0, ind = 0; j < 5; ++j) {
52  cov(j, j) = covariance(i)(ind++);
53  for (auto k = j + 1; k < 5; ++k)
54  cov(k, j) = cov(j, k) = covariance(i)(ind++);
55  }
56  }
57 };
58 
59 #endif // CUDADataFormats_Track_TrajectoryStateSOAT_H
Eigen::Matrix< double, 5, 5 > Matrix5d
eigenSoA::MatrixSoA< Vector15f, S > covariance
#define __host__
__host__ __device__ void copyToDense(V5 &v, M5 &cov, int32_t i) const
__host__ __device__ void copyFromCircle(V3 const &cp, M3 const &ccov, V2 const &lp, M2 const &lcov, float b, int32_t i)
Eigen::Matrix< float, 15, 1 > Vector15f
eigenSoA::MatrixSoA< Vector5f, S > state
static constexpr int32_t stride()
Eigen::Matrix< double, 5, 1 > Vector5d
__host__ __device__ void copyFromDense(V5 const &v, M5 const &cov, int32_t i)
double b
Definition: hdecay.h:118
#define __device__
Eigen::Matrix< float, 5, 1 > Vector5f