Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:05:23

0001 #include <type_traits>
0002 
0003 #include <alpaka/alpaka.hpp>
0004 
0005 #include "DataFormats/TrackSoA/interface/TrackDefinitions.h"
0006 #include "DataFormats/TrackSoA/interface/TracksDevice.h"
0007 #include "DataFormats/TrackSoA/interface/TracksHost.h"
0008 #include "DataFormats/TrackSoA/interface/alpaka/TracksSoACollection.h"
0009 #include "Geometry/CommonTopologies/interface/SimplePixelTopology.h"
0010 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0011 
0012 #include "TrackSoAHeterogeneous_test.h"
0013 
0014 using namespace reco;
0015 
0016 using Quality = pixelTrack::Quality;
0017 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0018   using namespace cms::alpakatools;
0019   namespace testTrackSoA {
0020 
0021     // Kernel which fills the TrackSoAView with data
0022     // to test writing to it
0023     template <typename TrackerTraits>
0024     class TestFillKernel {
0025     public:
0026       template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
0027       ALPAKA_FN_ACC void operator()(TAcc const& acc, TrackSoAView<TrackerTraits> tracks_view, int32_t nTracks) const {
0028         if (cms::alpakatools::once_per_grid(acc)) {
0029           tracks_view.nTracks() = nTracks;
0030         }
0031 
0032         for (int32_t j : uniform_elements(acc, nTracks)) {
0033           tracks_view[j].pt() = (float)j;
0034           tracks_view[j].eta() = (float)j;
0035           tracks_view[j].chi2() = (float)j;
0036           tracks_view[j].quality() = (Quality)(j % 256);
0037           tracks_view[j].nLayers() = j % 128;
0038           tracks_view.hitIndices().off[j] = j;
0039         }
0040       }
0041     };
0042 
0043     // Kernel which reads from the TrackSoAView to verify
0044     // that it was written correctly from the fill kernel
0045     template <typename TrackerTraits>
0046     class TestVerifyKernel {
0047     public:
0048       template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
0049       ALPAKA_FN_ACC void operator()(TAcc const& acc,
0050                                     TrackSoAConstView<TrackerTraits> tracks_view,
0051                                     int32_t nTracks) const {
0052         if (cms::alpakatools::once_per_grid(acc)) {
0053           ALPAKA_ASSERT(tracks_view.nTracks() == nTracks);
0054         }
0055         for (int32_t j : uniform_elements(acc, tracks_view.nTracks())) {
0056           ALPAKA_ASSERT(abs(tracks_view[j].pt() - (float)j) < .0001);
0057           ALPAKA_ASSERT(abs(tracks_view[j].eta() - (float)j) < .0001);
0058           ALPAKA_ASSERT(abs(tracks_view[j].chi2() - (float)j) < .0001);
0059           ALPAKA_ASSERT(tracks_view[j].quality() == (Quality)(j % 256));
0060           ALPAKA_ASSERT(tracks_view[j].nLayers() == j % 128);
0061           ALPAKA_ASSERT(tracks_view.hitIndices().off[j] == uint32_t(j));
0062         }
0063       }
0064     };
0065 
0066     // Host function which invokes the two kernels above
0067     template <typename TrackerTraits>
0068     void runKernels(TrackSoAView<TrackerTraits> tracks_view, Queue& queue) {
0069       int32_t tracks = 420;
0070       uint32_t items = 64;
0071       uint32_t groups = divide_up_by(tracks, items);
0072       auto workDiv = make_workdiv<Acc1D>(groups, items);
0073       alpaka::exec<Acc1D>(queue, workDiv, TestFillKernel<TrackerTraits>{}, tracks_view, tracks);
0074       alpaka::exec<Acc1D>(queue, workDiv, TestVerifyKernel<TrackerTraits>{}, tracks_view, tracks);
0075     }
0076 
0077     template void runKernels<pixelTopology::Phase1>(TrackSoAView<pixelTopology::Phase1> tracks_view, Queue& queue);
0078     template void runKernels<pixelTopology::Phase2>(TrackSoAView<pixelTopology::Phase2> tracks_view, Queue& queue);
0079 
0080   }  // namespace testTrackSoA
0081 }  // namespace ALPAKA_ACCELERATOR_NAMESPACE