Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include <cstdlib>
0002 #include <unistd.h>
0003 
0004 #include <alpaka/alpaka.hpp>
0005 
0006 #include "DataFormats/TrackingRecHitSoA/interface/TrackingRecHitsDevice.h"
0007 #include "DataFormats/TrackingRecHitSoA/interface/TrackingRecHitsHost.h"
0008 #include "DataFormats/TrackingRecHitSoA/interface/alpaka/TrackingRecHitsSoACollection.h"
0009 #include "FWCore/Utilities/interface/stringize.h"
0010 #include "Geometry/CommonTopologies/interface/SimplePixelTopology.h"
0011 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0012 #include "HeterogeneousCore/AlpakaInterface/interface/devices.h"
0013 #include "HeterogeneousCore/AlpakaInterface/interface/memory.h"
0014 #include "HeterogeneousCore/AlpakaInterface/interface/workdivision.h"
0015 
0016 #include "Hits_test.h"
0017 
0018 using namespace ALPAKA_ACCELERATOR_NAMESPACE;
0019 
0020 int main() {
0021   // Get the list of devices on the current platform
0022   auto const& devices = cms::alpakatools::devices<Platform>();
0023   if (devices.empty()) {
0024     std::cerr << "No devices available for the " EDM_STRINGIZE(ALPAKA_ACCELERATOR_NAMESPACE) " backend, "
0025       "the test will be skipped.\n";
0026     exit(EXIT_FAILURE);
0027   }
0028 
0029   // Run the test on each device
0030   for (const auto& device : devices) {
0031     Queue queue(device);
0032 
0033     // inner scope to deallocate memory before destroying the queue
0034     {
0035       uint32_t nHits = 2000;
0036       int32_t offset = 100;
0037       auto moduleStartH =
0038           cms::alpakatools::make_host_buffer<uint32_t[]>(queue, pixelTopology::Phase1::numberOfModules + 1);
0039       for (size_t i = 0; i < pixelTopology::Phase1::numberOfModules + 1; ++i) {
0040         moduleStartH[i] = i * 2;
0041       }
0042       auto moduleStartD =
0043           cms::alpakatools::make_device_buffer<uint32_t[]>(queue, pixelTopology::Phase1::numberOfModules + 1);
0044       alpaka::memcpy(queue, moduleStartD, moduleStartH);
0045       TrackingRecHitsSoACollection<pixelTopology::Phase1> tkhit(queue, nHits, offset, moduleStartD.data());
0046 
0047       testTrackingRecHitSoA::runKernels<pixelTopology::Phase1>(tkhit.view(), queue);
0048       tkhit.updateFromDevice(queue);
0049 
0050 #if defined ALPAKA_ACC_CPU_B_SEQ_T_SEQ_ENABLED or defined ALPAKA_ACC_CPU_B_TBB_T_SEQ_ENABLED
0051       // requires c++23 to make cms::alpakatools::CopyToHost compile using if constexpr
0052       // see https://www.open-std.org/jtc1/sc22/wg21/docs/papers/2022/p2593r0.html
0053       TrackingRecHitHost<pixelTopology::Phase1> const& host_collection = tkhit;
0054 #else
0055       TrackingRecHitHost<pixelTopology::Phase1> host_collection =
0056           cms::alpakatools::CopyToHost<TrackingRecHitDevice<pixelTopology::Phase1, Device> >::copyAsync(queue, tkhit);
0057 #endif
0058       // wait for the kernel and the potential copy to complete
0059       alpaka::wait(queue);
0060       assert(tkhit.nHits() == nHits);
0061       assert(tkhit.offsetBPIX2() == 22);  // set in the kernel
0062       assert(tkhit.nHits() == host_collection.nHits());
0063       assert(tkhit.offsetBPIX2() == host_collection.offsetBPIX2());
0064     }
0065   }
0066 
0067   return EXIT_SUCCESS;
0068 }