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/SiPixelDigiSoA/interface/SiPixelDigiErrorsDevice.h"
0004 #include "DataFormats/SiPixelDigiSoA/interface/SiPixelDigiErrorsHost.h"
0005 #include "DataFormats/SiPixelDigiSoA/interface/alpaka/SiPixelDigiErrorsSoACollection.h"
0006 #include "DataFormats/SiPixelRawData/interface/SiPixelErrorCompact.h"
0007 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0008 #include "HeterogeneousCore/AlpakaInterface/interface/workdivision.h"
0009 
0010 #include "DigiErrors_test.h"
0011 
0012 using namespace alpaka;
0013 
0014 namespace ALPAKA_ACCELERATOR_NAMESPACE::testDigisSoA {
0015 
0016   class TestFillKernel {
0017   public:
0018     ALPAKA_FN_ACC void operator()(Acc1D const& acc, SiPixelDigiErrorsSoAView digiErrors_view) const {
0019       for (uint32_t j : cms::alpakatools::uniform_elements(acc, digiErrors_view.metadata().size())) {
0020         digiErrors_view[j].pixelErrors().rawId = j;
0021         digiErrors_view[j].pixelErrors().word = j;
0022         digiErrors_view[j].pixelErrors().errorType = j;
0023         digiErrors_view[j].pixelErrors().fedId = j;
0024       }
0025     }
0026   };
0027 
0028   class TestVerifyKernel {
0029   public:
0030     ALPAKA_FN_ACC void operator()(Acc1D const& acc, SiPixelDigiErrorsSoAConstView digiErrors_view) const {
0031       for (uint32_t j : cms::alpakatools::uniform_elements(acc, digiErrors_view.metadata().size())) {
0032         ALPAKA_ASSERT_ACC(digiErrors_view[j].pixelErrors().rawId == j);
0033         ALPAKA_ASSERT_ACC(digiErrors_view[j].pixelErrors().word == j);
0034         ALPAKA_ASSERT_ACC(digiErrors_view[j].pixelErrors().errorType == j % 256);
0035         ALPAKA_ASSERT_ACC(digiErrors_view[j].pixelErrors().fedId == j % 256);
0036       }
0037     }
0038   };
0039 
0040   void runKernels(SiPixelDigiErrorsSoAView digiErrors_view, Queue& queue) {
0041     uint32_t items = 64;
0042     uint32_t groups = cms::alpakatools::divide_up_by(digiErrors_view.metadata().size(), items);
0043     auto workDiv = cms::alpakatools::make_workdiv<Acc1D>(groups, items);
0044     alpaka::exec<Acc1D>(queue, workDiv, TestFillKernel{}, digiErrors_view);
0045     alpaka::exec<Acc1D>(queue, workDiv, TestVerifyKernel{}, digiErrors_view);
0046   }
0047 
0048 }  // namespace ALPAKA_ACCELERATOR_NAMESPACE::testDigisSoA