Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-01-22 07:33:56

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       ALPAKA_FN_ACC void operator()(Acc1D const& acc, TrackSoAView<TrackerTraits> tracks_view, int32_t nTracks) const {
0027         if (cms::alpakatools::once_per_grid(acc)) {
0028           tracks_view.nTracks() = nTracks;
0029         }
0030 
0031         for (int32_t j : uniform_elements(acc, nTracks)) {
0032           tracks_view[j].pt() = (float)j;
0033           tracks_view[j].eta() = (float)j;
0034           tracks_view[j].chi2() = (float)j;
0035           tracks_view[j].quality() = (Quality)(j % 256);
0036           tracks_view[j].nLayers() = j % 128;
0037           tracks_view.hitIndices().off[j] = j;
0038         }
0039       }
0040     };
0041 
0042     // Kernel which reads from the TrackSoAView to verify
0043     // that it was written correctly from the fill kernel
0044     template <typename TrackerTraits>
0045     class TestVerifyKernel {
0046     public:
0047       ALPAKA_FN_ACC void operator()(Acc1D const& acc,
0048                                     TrackSoAConstView<TrackerTraits> tracks_view,
0049                                     int32_t nTracks) const {
0050         if (cms::alpakatools::once_per_grid(acc)) {
0051           ALPAKA_ASSERT(tracks_view.nTracks() == nTracks);
0052         }
0053         for (int32_t j : uniform_elements(acc, tracks_view.nTracks())) {
0054           ALPAKA_ASSERT(abs(tracks_view[j].pt() - (float)j) < .0001);
0055           ALPAKA_ASSERT(abs(tracks_view[j].eta() - (float)j) < .0001);
0056           ALPAKA_ASSERT(abs(tracks_view[j].chi2() - (float)j) < .0001);
0057           ALPAKA_ASSERT(tracks_view[j].quality() == (Quality)(j % 256));
0058           ALPAKA_ASSERT(tracks_view[j].nLayers() == j % 128);
0059           ALPAKA_ASSERT(tracks_view.hitIndices().off[j] == uint32_t(j));
0060         }
0061       }
0062     };
0063 
0064     // Host function which invokes the two kernels above
0065     template <typename TrackerTraits>
0066     void runKernels(TrackSoAView<TrackerTraits> tracks_view, Queue& queue) {
0067       int32_t tracks = 420;
0068       uint32_t items = 64;
0069       uint32_t groups = divide_up_by(tracks, items);
0070       auto workDiv = make_workdiv<Acc1D>(groups, items);
0071       alpaka::exec<Acc1D>(queue, workDiv, TestFillKernel<TrackerTraits>{}, tracks_view, tracks);
0072       alpaka::exec<Acc1D>(queue, workDiv, TestVerifyKernel<TrackerTraits>{}, tracks_view, tracks);
0073     }
0074 
0075     template void runKernels<pixelTopology::Phase1>(TrackSoAView<pixelTopology::Phase1> tracks_view, Queue& queue);
0076     template void runKernels<pixelTopology::Phase2>(TrackSoAView<pixelTopology::Phase2> tracks_view, Queue& queue);
0077 
0078   }  // namespace testTrackSoA
0079 }  // namespace ALPAKA_ACCELERATOR_NAMESPACE