File indexing completed on 2025-07-03 00:42:11
0001 #include <Eigen/Dense>
0002
0003 #include "DataFormats/TrackSoA/interface/TracksSoA.h"
0004 #include "DataFormats/TrackSoA/interface/alpaka/TrackUtilities.h"
0005 #include "Geometry/CommonTopologies/interface/SimplePixelTopology.h"
0006 #include "HeterogeneousCore/AlpakaInterface/interface/workdivision.h"
0007
0008 #include "TrajectoryStateSoA_t.h"
0009
0010 using Vector5d = Eigen::Matrix<double, 5, 1>;
0011 using Matrix5d = Eigen::Matrix<double, 5, 5>;
0012
0013 using namespace cms::alpakatools;
0014
0015 namespace ALPAKA_ACCELERATOR_NAMESPACE::test {
0016
0017 namespace {
0018
0019 ALPAKA_FN_ACC Matrix5d buildCovariance(Vector5d const& e) {
0020 Matrix5d cov;
0021 for (int i = 0; i < 5; ++i)
0022 cov(i, i) = e(i) * e(i);
0023 for (int i = 0; i < 5; ++i) {
0024 for (int j = 0; j < i; ++j) {
0025
0026 double v = 0.3 * std::sqrt(cov(i, i) * cov(j, j));
0027 cov(i, j) = (i + j) % 2 ? -0.4 * v : 0.1 * v;
0028 cov(j, i) = cov(i, j);
0029 }
0030 }
0031 return cov;
0032 }
0033
0034 template <typename TrackerTraits>
0035 struct TestTrackSoA {
0036 ALPAKA_FN_ACC void operator()(Acc1D const& acc, reco::TrackSoAView tracks) const {
0037 Vector5d par0;
0038 par0 << 0.2, 0.1, 3.5, 0.8, 0.1;
0039 Vector5d e0;
0040 e0 << 0.01, 0.01, 0.035, -0.03, -0.01;
0041 Matrix5d cov0 = buildCovariance(e0);
0042
0043 for (auto i : uniform_elements(acc, tracks.metadata().size())) {
0044 reco::copyFromDense(tracks, par0, cov0, i);
0045 Vector5d par1;
0046 Matrix5d cov1;
0047 reco::copyToDense(tracks, par1, cov1, i);
0048 Vector5d deltaV = par1 - par0;
0049 Matrix5d deltaM = cov1 - cov0;
0050 for (int j = 0; j < 5; ++j) {
0051 ALPAKA_ASSERT(std::abs(deltaV(j)) < 1.e-5);
0052 for (int k = j; k < 5; ++k) {
0053 ALPAKA_ASSERT(cov0(k, j) == cov0(j, k));
0054 ALPAKA_ASSERT(cov1(k, j) == cov1(j, k));
0055 ALPAKA_ASSERT(std::abs(deltaM(k, j)) < 1.e-5);
0056 }
0057 }
0058 }
0059 }
0060 };
0061
0062 }
0063
0064 template <typename TrackerTraits>
0065 void testTrackSoA(Queue& queue, ::reco::TrackSoAView& tracks) {
0066 auto grid = make_workdiv<Acc1D>(1, 64);
0067 alpaka::exec<Acc1D>(queue, grid, TestTrackSoA<TrackerTraits>{}, tracks);
0068 }
0069
0070 template void testTrackSoA<pixelTopology::Phase1>(Queue& queue, reco::TrackSoAView& tracks);
0071 template void testTrackSoA<pixelTopology::Phase2>(Queue& queue, reco::TrackSoAView& tracks);
0072
0073 }