Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include <alpaka/alpaka.hpp>
0002 
0003 #include "DataFormats/SiPixelClusterSoA/interface/SiPixelClustersDevice.h"
0004 #include "DataFormats/SiPixelClusterSoA/interface/SiPixelClustersHost.h"
0005 #include "DataFormats/SiPixelClusterSoA/interface/alpaka/SiPixelClustersSoACollection.h"
0006 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0007 #include "HeterogeneousCore/AlpakaInterface/interface/workdivision.h"
0008 
0009 #include "Clusters_test.h"
0010 
0011 using namespace alpaka;
0012 
0013 namespace ALPAKA_ACCELERATOR_NAMESPACE::testClusterSoA {
0014 
0015   class TestFillKernel {
0016   public:
0017     ALPAKA_FN_ACC void operator()(Acc1D const& acc, SiPixelClustersSoAView clust_view) const {
0018       for (int32_t j : cms::alpakatools::uniform_elements(acc, clust_view.metadata().size())) {
0019         clust_view[j].moduleStart() = j;
0020         clust_view[j].clusInModule() = j * 2;
0021         clust_view[j].moduleId() = j * 3;
0022         clust_view[j].clusModuleStart() = j * 4;
0023       }
0024     }
0025   };
0026 
0027   class TestVerifyKernel {
0028   public:
0029     ALPAKA_FN_ACC void operator()(Acc1D const& acc, SiPixelClustersSoAConstView clust_view) const {
0030       for (uint32_t j : cms::alpakatools::uniform_elements(acc, clust_view.metadata().size())) {
0031         ALPAKA_ASSERT_ACC(clust_view[j].moduleStart() == j);
0032         ALPAKA_ASSERT_ACC(clust_view[j].clusInModule() == j * 2);
0033         ALPAKA_ASSERT_ACC(clust_view[j].moduleId() == j * 3);
0034         ALPAKA_ASSERT_ACC(clust_view[j].clusModuleStart() == j * 4);
0035       }
0036     }
0037   };
0038 
0039   void runKernels(SiPixelClustersSoAView clust_view, Queue& queue) {
0040     uint32_t items = 64;
0041     uint32_t groups = cms::alpakatools::divide_up_by(clust_view.metadata().size(), items);
0042     auto workDiv = cms::alpakatools::make_workdiv<Acc1D>(groups, items);
0043     alpaka::exec<Acc1D>(queue, workDiv, TestFillKernel{}, clust_view);
0044     alpaka::exec<Acc1D>(queue, workDiv, TestVerifyKernel{}, clust_view);
0045   }
0046 
0047 }  // namespace ALPAKA_ACCELERATOR_NAMESPACE::testClusterSoA