Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-20 02:31:45

0001 #include <type_traits>
0002 
0003 #include <alpaka/alpaka.hpp>
0004 
0005 #include "DataFormats/TrackingRecHitSoA/interface/TrackingRecHitsDevice.h"
0006 #include "DataFormats/TrackingRecHitSoA/interface/TrackingRecHitsSoA.h"
0007 #include "DataFormats/TrackingRecHitSoA/interface/alpaka/TrackingRecHitsSoACollection.h"
0008 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0009 #include "HeterogeneousCore/AlpakaInterface/interface/workdivision.h"
0010 
0011 #include "Hits_test.h"
0012 
0013 using namespace alpaka;
0014 
0015 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0016 
0017   using namespace cms::alpakatools;
0018   namespace testTrackingRecHitSoA {
0019 
0020     template <typename TrackerTraits>
0021     struct TestFillKernel {
0022       template <typename TAcc, typename = std::enable_if_t<isAccelerator<TAcc>>>
0023       ALPAKA_FN_ACC void operator()(TAcc const& acc, TrackingRecHitSoAView<TrackerTraits> soa) const {
0024         const uint32_t i(alpaka::getIdx<alpaka::Grid, alpaka::Blocks>(acc)[0u]);
0025         const uint32_t j(alpaka::getIdx<alpaka::Block, alpaka::Threads>(acc)[0u]);
0026 
0027         if (cms::alpakatools::once_per_grid(acc)) {
0028           soa.offsetBPIX2() = 22;
0029           soa[10].xLocal() = 1.11;
0030         }
0031 
0032         soa[i].iphi() = i % 10;
0033         soa.hitsLayerStart()[j] = j;
0034       }
0035     };
0036 
0037     template <typename TrackerTraits>
0038     struct ShowKernel {
0039       template <typename TAcc, typename = std::enable_if_t<isAccelerator<TAcc>>>
0040       ALPAKA_FN_ACC void operator()(TAcc const& acc, TrackingRecHitSoAConstView<TrackerTraits> soa) const {
0041         if (cms::alpakatools::once_per_grid(acc)) {
0042           printf("nbins = %d\n", soa.phiBinner().nbins());
0043           printf("offsetBPIX = %d\n", soa.offsetBPIX2());
0044           printf("nHits = %d\n", soa.metadata().size());
0045           //printf("hitsModuleStart[28] = %d\n", soa[28].hitsModuleStart());
0046         }
0047 
0048         // can be increased to soa.nHits() for debugging
0049         for (uint32_t i : cms::alpakatools::uniform_elements(acc, 10)) {
0050           printf("iPhi %d -> %d\n", i, soa[i].iphi());
0051         }
0052       }
0053     };
0054 
0055     template <typename TrackerTraits>
0056     void runKernels(TrackingRecHitSoAView<TrackerTraits>& view, Queue& queue) {
0057       uint32_t items = 64;
0058       uint32_t groups = divide_up_by(view.metadata().size(), items);
0059       auto workDiv = make_workdiv<Acc1D>(groups, items);
0060       alpaka::exec<Acc1D>(queue, workDiv, TestFillKernel<TrackerTraits>{}, view);
0061       alpaka::exec<Acc1D>(queue, workDiv, ShowKernel<TrackerTraits>{}, view);
0062     }
0063 
0064     template void runKernels<pixelTopology::Phase1>(TrackingRecHitSoAView<pixelTopology::Phase1>& view, Queue& queue);
0065     template void runKernels<pixelTopology::Phase2>(TrackingRecHitSoAView<pixelTopology::Phase2>& view, Queue& queue);
0066 
0067   }  // namespace testTrackingRecHitSoA
0068 }  // namespace ALPAKA_ACCELERATOR_NAMESPACE