Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-04 04:35:00

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           // this makes the matrix positive defined
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       using Utils = TracksUtilities<TrackerTraits>;
0037 
0038       ALPAKA_FN_ACC void operator()(Acc1D const& acc, reco::TrackSoAView<TrackerTraits> tracks) const {
0039         Vector5d par0;
0040         par0 << 0.2, 0.1, 3.5, 0.8, 0.1;
0041         Vector5d e0;
0042         e0 << 0.01, 0.01, 0.035, -0.03, -0.01;
0043         Matrix5d cov0 = buildCovariance(e0);
0044 
0045         for (auto i : uniform_elements(acc, tracks.metadata().size())) {
0046           Utils::copyFromDense(tracks, par0, cov0, i);
0047           Vector5d par1;
0048           Matrix5d cov1;
0049           Utils::copyToDense(tracks, par1, cov1, i);
0050           Vector5d deltaV = par1 - par0;
0051           Matrix5d deltaM = cov1 - cov0;
0052           for (int j = 0; j < 5; ++j) {
0053             ALPAKA_ASSERT(std::abs(deltaV(j)) < 1.e-5);
0054             for (int k = j; k < 5; ++k) {
0055               ALPAKA_ASSERT(cov0(k, j) == cov0(j, k));
0056               ALPAKA_ASSERT(cov1(k, j) == cov1(j, k));
0057               ALPAKA_ASSERT(std::abs(deltaM(k, j)) < 1.e-5);
0058             }
0059           }
0060         }
0061       }
0062     };
0063 
0064   }  // namespace
0065 
0066   template <typename TrackerTraits>
0067   void testTrackSoA(Queue& queue, reco::TrackSoAView<TrackerTraits>& tracks) {
0068     auto grid = make_workdiv<Acc1D>(1, 64);
0069     alpaka::exec<Acc1D>(queue, grid, TestTrackSoA<TrackerTraits>{}, tracks);
0070   }
0071 
0072   template void testTrackSoA<pixelTopology::Phase1>(Queue& queue, reco::TrackSoAView<pixelTopology::Phase1>& tracks);
0073   template void testTrackSoA<pixelTopology::Phase2>(Queue& queue, reco::TrackSoAView<pixelTopology::Phase2>& tracks);
0074 
0075 }  // namespace ALPAKA_ACCELERATOR_NAMESPACE::test