Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-07-03 00:42:11

0001 #include <cstdlib>
0002 #include <unistd.h>
0003 
0004 #include <alpaka/alpaka.hpp>
0005 
0006 #include "DataFormats/SiPixelClusterSoA/interface/alpaka/SiPixelClustersSoACollection.h"
0007 #include "DataFormats/TrackingRecHitSoA/interface/TrackingRecHitsDevice.h"
0008 #include "DataFormats/TrackingRecHitSoA/interface/TrackingRecHitsHost.h"
0009 #include "DataFormats/TrackingRecHitSoA/interface/alpaka/TrackingRecHitsSoACollection.h"
0010 #include "FWCore/Utilities/interface/stringize.h"
0011 #include "Geometry/CommonTopologies/interface/SimplePixelTopology.h"
0012 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0013 #include "HeterogeneousCore/AlpakaInterface/interface/devices.h"
0014 #include "HeterogeneousCore/AlpakaInterface/interface/memory.h"
0015 #include "HeterogeneousCore/AlpakaInterface/interface/workdivision.h"
0016 
0017 #include "Hits_test.h"
0018 
0019 using namespace ALPAKA_ACCELERATOR_NAMESPACE;
0020 using namespace ALPAKA_ACCELERATOR_NAMESPACE::reco;
0021 
0022 int main() {
0023   // Get the list of devices on the current platform
0024   auto const& devices = cms::alpakatools::devices<Platform>();
0025   if (devices.empty()) {
0026     std::cerr << "No devices available for the " EDM_STRINGIZE(ALPAKA_ACCELERATOR_NAMESPACE) " backend, "
0027       "the test will be skipped.\n";
0028     exit(EXIT_FAILURE);
0029   }
0030 
0031   // Run the test on each device
0032   for (const auto& device : devices) {
0033     Queue queue(device);
0034 
0035     // inner scope to deallocate memory before destroying the queue
0036     {
0037       uint32_t nHits = 2000;
0038       int32_t offset = 100;
0039       uint32_t nModules = 200;
0040 
0041       SiPixelClustersSoACollection clusters(nModules, queue);
0042       clusters.setNClusters(nHits, offset);
0043 
0044       auto moduleStartH = cms::alpakatools::make_host_buffer<uint32_t[]>(queue, nModules + 1);
0045 
0046       for (size_t i = 0; i < nModules + 1; ++i) {
0047         moduleStartH[i] = i * 2;
0048       }
0049 
0050       auto hitsX = cms::alpakatools::make_host_buffer<float[]>(queue, nHits);
0051       for (size_t i = 0; i < nHits; ++i) {
0052         hitsX[i] = float(i) * 2;
0053       }
0054 
0055       auto moduleStartD = cms::alpakatools::make_device_view<uint32_t>(queue, clusters.view().clusModuleStart(), nHits);
0056       alpaka::memcpy(queue, moduleStartD, moduleStartH);
0057 
0058       TrackingRecHitsSoACollection tkhit(queue, clusters);
0059 
0060       // exercise the copy of a full column (on device)
0061       auto hitXD = cms::alpakatools::make_device_view<float>(queue, tkhit.view().xLocal(), nHits);
0062       alpaka::memcpy(queue, hitXD, hitsX);
0063 
0064       // exercise the memset of a colum (on device)
0065       auto hitYD = cms::alpakatools::make_device_view<float>(queue, tkhit.view().yGlobal(), nHits);
0066       constexpr float constYG = -14.0458;
0067       std::vector<float> constYV(nHits, constYG);
0068       auto constYGV_v = cms::alpakatools::make_host_view<float>(constYV.data(), nHits);
0069       alpaka::memcpy(queue, hitYD, constYGV_v);
0070 
0071       testTrackingRecHitSoA::runKernels(tkhit.view(), tkhit.view<::reco::HitModuleSoA>(), queue);
0072       tkhit.updateFromDevice(queue);
0073 
0074 #if defined ALPAKA_ACC_CPU_B_SEQ_T_SEQ_ENABLED or defined ALPAKA_ACC_CPU_B_TBB_T_SEQ_ENABLED
0075       // requires c++23 to make cms::alpakatools::CopyToHost compile using if constexpr
0076       // see https://www.open-std.org/jtc1/sc22/wg21/docs/papers/2022/p2593r0.html
0077       ::reco::TrackingRecHitHost const& host_collection = tkhit;
0078 #else
0079       ::reco::TrackingRecHitHost host_collection =
0080           cms::alpakatools::CopyToHost<::reco::TrackingRecHitDevice<Device>>::copyAsync(queue, tkhit);
0081       alpaka::wait(queue);
0082 #endif
0083 
0084       alpaka::QueueCpuBlocking queue_host{cms::alpakatools::host()};
0085 
0086       ::reco::TrackingRecHitHost host_collection_2(cms::alpakatools::host(), nHits, nModules);
0087 
0088       // exercise the memset of a colum (on host)
0089       auto hitLYH = cms::alpakatools::make_host_view<float>(host_collection_2.view().yLocal(), nHits);
0090       constexpr float constYL = -27.0855;
0091       std::vector<float> constYLV(nHits, constYL);
0092       auto constYL_v = cms::alpakatools::make_host_view<float>(constYLV.data(), nHits);
0093       alpaka::memcpy(queue_host, hitLYH, constYL_v);
0094 
0095       // wait for the copy above to complete
0096       alpaka::wait(queue_host);
0097 
0098       assert(host_collection.view().xLocal()[12] == 24.);
0099       assert(host_collection.view().yGlobal()[int(nHits / 2)] == constYG);
0100       assert(host_collection_2.view().yLocal()[nHits - 1] == constYL);
0101 
0102       assert(tkhit.nHits() == nHits);
0103       assert(tkhit.offsetBPIX2() == 22);  // set in the kernel
0104       assert(tkhit.nHits() == host_collection.nHits());
0105       assert(tkhit.offsetBPIX2() == host_collection.offsetBPIX2());
0106     }
0107   }
0108 
0109   return EXIT_SUCCESS;
0110 }