Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-06-26 23:26:38

0001 #include <Eigen/Core>
0002 #include <Eigen/Dense>
0003 
0004 #include <alpaka/alpaka.hpp>
0005 
0006 #define CATCH_CONFIG_MAIN
0007 #include <catch.hpp>
0008 
0009 #include "DataFormats/SoATemplate/interface/SoALayout.h"
0010 #include "DataFormats/Portable/interface/PortableCollection.h"
0011 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0012 #include "HeterogeneousCore/AlpakaInterface/interface/memory.h"
0013 #include "HeterogeneousCore/AlpakaInterface/interface/workdivision.h"
0014 
0015 using namespace ALPAKA_ACCELERATOR_NAMESPACE;
0016 
0017 GENERATE_SOA_LAYOUT(SoAPositionTemplate,
0018                     SOA_COLUMN(float, x),
0019                     SOA_COLUMN(float, y),
0020                     SOA_COLUMN(float, z),
0021                     SOA_SCALAR(int, detectorType))
0022 
0023 using SoAPosition = SoAPositionTemplate<>;
0024 using SoAPositionView = SoAPosition::View;
0025 using SoAPositionConstView = SoAPosition::ConstView;
0026 
0027 GENERATE_SOA_LAYOUT(SoAPCATemplate,
0028                     SOA_COLUMN(float, eigenvalues),
0029                     SOA_COLUMN(float, eigenvector_1),
0030                     SOA_COLUMN(float, eigenvector_2),
0031                     SOA_COLUMN(float, eigenvector_3),
0032                     SOA_EIGEN_COLUMN(Eigen::Vector3d, candidateDirection))
0033 
0034 using SoAPCA = SoAPCATemplate<>;
0035 using SoAPCAView = SoAPCA::View;
0036 using SoAPCAConstView = SoAPCA::ConstView;
0037 
0038 GENERATE_SOA_LAYOUT(GenericSoATemplate,
0039                     SOA_COLUMN(float, x),
0040                     SOA_COLUMN(float, y),
0041                     SOA_COLUMN(float, z),
0042                     SOA_EIGEN_COLUMN(Eigen::Vector3d, candidateDirection))
0043 
0044 using GenericSoA = GenericSoATemplate<cms::soa::CacheLineSize::IntelCPU>;
0045 using GenericSoAView = GenericSoA::View;
0046 using GenericSoAConstView = GenericSoA::ConstView;
0047 
0048 // Kernel for filling the SoA
0049 struct FillSoA {
0050   template <typename TAcc, typename PositionView, typename PCAView>
0051   ALPAKA_FN_ACC void operator()(TAcc const& acc, PositionView positionView, PCAView pcaView) const {
0052     constexpr float interval = 0.01f;
0053     if (cms::alpakatools::once_per_grid(acc))
0054       positionView.detectorType() = 1;
0055 
0056     for (auto local_idx : cms::alpakatools::uniform_elements(acc, positionView.metadata().size())) {
0057       positionView[local_idx].x() = static_cast<float>(local_idx);
0058       positionView[local_idx].y() = static_cast<float>(local_idx) * 2.0f;
0059       positionView[local_idx].z() = static_cast<float>(local_idx) * 3.0f;
0060 
0061       pcaView[local_idx].eigenvector_1() = positionView[local_idx].x() / interval;
0062       pcaView[local_idx].eigenvector_2() = positionView[local_idx].y() / interval;
0063       pcaView[local_idx].eigenvector_3() = positionView[local_idx].z() / interval;
0064       pcaView[local_idx].candidateDirection()(0) = positionView[local_idx].x() / interval;
0065       pcaView[local_idx].candidateDirection()(1) = positionView[local_idx].y() / interval;
0066       pcaView[local_idx].candidateDirection()(2) = positionView[local_idx].z() / interval;
0067     }
0068   }
0069 };
0070 
0071 TEST_CASE("Deep copy from SoA Generic View") {
0072   auto const& devices = cms::alpakatools::devices<Platform>();
0073   if (devices.empty()) {
0074     FAIL("No devices available for the " EDM_STRINGIZE(ALPAKA_ACCELERATOR_NAMESPACE) " backend, "
0075         "the test will be skipped.");
0076   }
0077 
0078   auto devHost = alpaka::getDevByIdx(alpaka::PlatformCpu{}, 0u);
0079 
0080   for (auto const& device : cms::alpakatools::devices<Platform>()) {
0081     std::cout << "Running on " << alpaka::getName(device) << std::endl;
0082 
0083     Queue queue(device);
0084 
0085     // common number of elements for the SoAs
0086     const std::size_t elems = 10;
0087 
0088     // Portable Collections
0089     PortableCollection<SoAPosition, Device> positionCollection(elems, queue);
0090     PortableCollection<SoAPCA, Device> pcaCollection(elems, queue);
0091 
0092     // Portable Collection Views
0093     SoAPositionView& positionCollectionView = positionCollection.view();
0094     SoAPCAView& pcaCollectionView = pcaCollection.view();
0095     // Portable Collection ConstViews
0096     const SoAPositionConstView& positionCollectionConstView = positionCollection.const_view();
0097     const SoAPCAConstView& pcaCollectionConstView = pcaCollection.const_view();
0098 
0099     // fill up
0100     auto blockSize = 64;
0101     auto numberOfBlocks = cms::alpakatools::divide_up_by(elems, blockSize);
0102 
0103     const auto workDiv = cms::alpakatools::make_workdiv<Acc1D>(numberOfBlocks, blockSize);
0104 
0105     alpaka::exec<Acc1D>(queue, workDiv, FillSoA{}, positionCollectionView, pcaCollectionView);
0106 
0107     alpaka::wait(queue);
0108 
0109     SECTION("Deep copy the View host to host and device to device") {
0110       // addresses and size of the SoA columns
0111       const auto posRecs = positionCollectionView.records();
0112       const auto pcaRecs = pcaCollectionView.records();
0113 
0114       // building the View with runtime check for the size
0115       GenericSoAView genericView(posRecs.x(), posRecs.y(), posRecs.z(), pcaRecs.candidateDirection());
0116 
0117       // Check for equality of memory addresses
0118       REQUIRE(genericView.metadata().addressOf_x() == positionCollectionView.metadata().addressOf_x());
0119       REQUIRE(genericView.metadata().addressOf_y() == positionCollectionView.metadata().addressOf_y());
0120       REQUIRE(genericView.metadata().addressOf_z() == positionCollectionView.metadata().addressOf_z());
0121       REQUIRE(genericView.metadata().addressOf_candidateDirection() ==
0122               pcaCollectionView.metadata().addressOf_candidateDirection());
0123 
0124       // PortableCollection that will host the aggregated columns
0125       PortableCollection<GenericSoA, Device> genericCollection(elems, queue);
0126       genericCollection.deepCopy(genericView, queue);
0127 
0128       // Check for inequality of memory addresses[`View` section](../../DataFormats/SoATemplate/README.md#view)
0129       REQUIRE(genericCollection.view().metadata().addressOf_x() != positionCollectionView.metadata().addressOf_x());
0130       REQUIRE(genericCollection.view().metadata().addressOf_y() != positionCollectionView.metadata().addressOf_y());
0131       REQUIRE(genericCollection.view().metadata().addressOf_z() != positionCollectionView.metadata().addressOf_z());
0132       REQUIRE(genericCollection.view().metadata().addressOf_candidateDirection() !=
0133               pcaCollectionView.metadata().addressOf_candidateDirection());
0134     }
0135 
0136     SECTION("Deep copy the ConstView host to host and device to device") {
0137       // addresses and size of the SoA columns
0138       const auto posRecs = positionCollectionConstView.records();
0139       const auto pcaRecs = pcaCollectionConstView.records();
0140 
0141       // building the View with runtime check for the size
0142       GenericSoAConstView genericConstView(posRecs.x(), posRecs.y(), posRecs.z(), pcaRecs.candidateDirection());
0143 
0144       // Check for equality of memory addresses
0145       REQUIRE(genericConstView.metadata().addressOf_x() == positionCollectionView.metadata().addressOf_x());
0146       REQUIRE(genericConstView.metadata().addressOf_y() == positionCollectionView.metadata().addressOf_y());
0147       REQUIRE(genericConstView.metadata().addressOf_z() == positionCollectionView.metadata().addressOf_z());
0148       REQUIRE(genericConstView.metadata().addressOf_candidateDirection() ==
0149               pcaCollectionView.metadata().addressOf_candidateDirection());
0150 
0151       // PortableCollection that will host the aggregated columns
0152       PortableCollection<GenericSoA, Device> genericCollection(elems, queue);
0153       genericCollection.deepCopy(genericConstView, queue);
0154 
0155       // Check for inequality of memory addresses
0156       REQUIRE(genericCollection.view().metadata().addressOf_x() != positionCollectionView.metadata().addressOf_x());
0157       REQUIRE(genericCollection.view().metadata().addressOf_y() != positionCollectionView.metadata().addressOf_y());
0158       REQUIRE(genericCollection.view().metadata().addressOf_z() != positionCollectionView.metadata().addressOf_z());
0159       REQUIRE(genericCollection.view().metadata().addressOf_candidateDirection() !=
0160               pcaCollectionView.metadata().addressOf_candidateDirection());
0161 
0162       // Check for correctness of the copy
0163       PortableHostCollection<GenericSoA> genericHostCollection(elems, queue);
0164       PortableHostCollection<SoAPosition> positionHostCollection(elems, queue);
0165       PortableHostCollection<SoAPCA> pcaHostCollection(elems, queue);
0166 
0167       alpaka::memcpy(queue, genericHostCollection.buffer(), genericCollection.buffer());
0168       alpaka::memcpy(queue, positionHostCollection.buffer(), positionCollection.buffer());
0169       alpaka::memcpy(queue, pcaHostCollection.buffer(), pcaCollection.buffer());
0170 
0171       alpaka::wait(queue);
0172 
0173       const GenericSoAConstView& genericViewHostCollection = genericHostCollection.const_view();
0174       const SoAPositionConstView& positionViewHostCollection = positionHostCollection.const_view();
0175       const SoAPCAConstView& pcaViewHostCollection = pcaHostCollection.const_view();
0176 
0177       for (size_t i = 0; i < elems; i++) {
0178         REQUIRE(genericViewHostCollection[i].x() == positionViewHostCollection[i].x());
0179         REQUIRE(genericViewHostCollection[i].y() == positionViewHostCollection[i].y());
0180         REQUIRE(genericViewHostCollection[i].z() == positionViewHostCollection[i].z());
0181         REQUIRE(genericViewHostCollection[i].candidateDirection()(0) ==
0182                 pcaViewHostCollection[i].candidateDirection()(0));
0183         REQUIRE(genericViewHostCollection[i].candidateDirection()(1) ==
0184                 pcaViewHostCollection[i].candidateDirection()(1));
0185         REQUIRE(genericViewHostCollection[i].candidateDirection()(2) ==
0186                 pcaViewHostCollection[i].candidateDirection()(2));
0187       }
0188     }
0189 
0190     SECTION("Deep copy the ConstView device to host") {
0191       // addresses and size of the SoA columns
0192       const auto posRecs = positionCollectionConstView.records();
0193       const auto pcaRecs = pcaCollectionConstView.records();
0194 
0195       // building the View with runtime check for the size
0196       GenericSoAConstView genericConstView(posRecs.x(), posRecs.y(), posRecs.z(), pcaRecs.candidateDirection());
0197 
0198       // Check for equality of memory addresses
0199       REQUIRE(genericConstView.metadata().addressOf_x() == positionCollectionView.metadata().addressOf_x());
0200       REQUIRE(genericConstView.metadata().addressOf_y() == positionCollectionView.metadata().addressOf_y());
0201       REQUIRE(genericConstView.metadata().addressOf_z() == positionCollectionView.metadata().addressOf_z());
0202       REQUIRE(genericConstView.metadata().addressOf_candidateDirection() ==
0203               pcaCollectionView.metadata().addressOf_candidateDirection());
0204 
0205       // PortableCollection that will host the aggregated columns
0206       PortableHostCollection<GenericSoA> genericCollection(elems, queue);
0207       genericCollection.deepCopy(genericConstView, queue);
0208 
0209       // Check for inequality of memory addresses
0210       REQUIRE(genericCollection.view().metadata().addressOf_x() != positionCollectionView.metadata().addressOf_x());
0211       REQUIRE(genericCollection.view().metadata().addressOf_y() != positionCollectionView.metadata().addressOf_y());
0212       REQUIRE(genericCollection.view().metadata().addressOf_z() != positionCollectionView.metadata().addressOf_z());
0213       REQUIRE(genericCollection.view().metadata().addressOf_candidateDirection() !=
0214               pcaCollectionView.metadata().addressOf_candidateDirection());
0215 
0216       // Check for correctness of the copy
0217       PortableHostCollection<SoAPosition> positionHostCollection(elems, queue);
0218       PortableHostCollection<SoAPCA> pcaHostCollection(elems, queue);
0219 
0220       alpaka::memcpy(queue, positionHostCollection.buffer(), positionCollection.buffer());
0221       alpaka::memcpy(queue, pcaHostCollection.buffer(), pcaCollection.buffer());
0222 
0223       alpaka::wait(queue);
0224 
0225       const GenericSoAConstView& genericViewCollection = genericCollection.const_view();
0226       const SoAPositionConstView& positionViewHostCollection = positionHostCollection.const_view();
0227       const SoAPCAConstView& pcaViewHostCollection = pcaHostCollection.const_view();
0228 
0229       for (size_t i = 0; i < elems; i++) {
0230         REQUIRE(genericViewCollection[i].x() == positionViewHostCollection[i].x());
0231         REQUIRE(genericViewCollection[i].y() == positionViewHostCollection[i].y());
0232         REQUIRE(genericViewCollection[i].z() == positionViewHostCollection[i].z());
0233         REQUIRE(genericViewCollection[i].candidateDirection()(0) == pcaViewHostCollection[i].candidateDirection()(0));
0234         REQUIRE(genericViewCollection[i].candidateDirection()(1) == pcaViewHostCollection[i].candidateDirection()(1));
0235         REQUIRE(genericViewCollection[i].candidateDirection()(2) == pcaViewHostCollection[i].candidateDirection()(2));
0236       }
0237     }
0238 
0239     SECTION("Deep copy the ConstView host to device") {
0240       PortableHostCollection<SoAPosition> positionHostCollection(elems, queue);
0241       PortableHostCollection<SoAPCA> pcaHostCollection(elems, queue);
0242 
0243       alpaka::memcpy(queue, positionHostCollection.buffer(), positionCollection.buffer());
0244       alpaka::memcpy(queue, pcaHostCollection.buffer(), pcaCollection.buffer());
0245 
0246       const SoAPositionConstView& positionViewHostCollection = positionHostCollection.const_view();
0247       const SoAPCAConstView& pcaViewHostCollection = pcaHostCollection.const_view();
0248 
0249       // addresses and size of the SoA columns
0250       const auto posRecs = positionViewHostCollection.records();
0251       const auto pcaRecs = pcaViewHostCollection.records();
0252 
0253       // building the View with runtime check for the size
0254       GenericSoAConstView genericConstView(posRecs.x(), posRecs.y(), posRecs.z(), pcaRecs.candidateDirection());
0255 
0256       // Check for equality of memory addresses
0257       REQUIRE(genericConstView.metadata().addressOf_x() == positionViewHostCollection.metadata().addressOf_x());
0258       REQUIRE(genericConstView.metadata().addressOf_y() == positionViewHostCollection.metadata().addressOf_y());
0259       REQUIRE(genericConstView.metadata().addressOf_z() == positionViewHostCollection.metadata().addressOf_z());
0260       REQUIRE(genericConstView.metadata().addressOf_candidateDirection() ==
0261               pcaViewHostCollection.metadata().addressOf_candidateDirection());
0262 
0263       // PortableCollection that will host the aggregated columns
0264       PortableCollection<GenericSoA, Device> genericCollection(elems, queue);
0265       genericCollection.deepCopy(genericConstView, queue);
0266 
0267       // Check for inequality of memory addresses
0268       REQUIRE(genericCollection.view().metadata().addressOf_x() != positionViewHostCollection.metadata().addressOf_x());
0269       REQUIRE(genericCollection.view().metadata().addressOf_y() != positionViewHostCollection.metadata().addressOf_y());
0270       REQUIRE(genericCollection.view().metadata().addressOf_z() != positionViewHostCollection.metadata().addressOf_z());
0271       REQUIRE(genericCollection.view().metadata().addressOf_candidateDirection() !=
0272               pcaViewHostCollection.metadata().addressOf_candidateDirection());
0273 
0274       // Check for correctness of the copy
0275       PortableHostCollection<GenericSoA> genericHostCollection(elems, queue);
0276 
0277       alpaka::memcpy(queue, genericHostCollection.buffer(), genericCollection.buffer());
0278 
0279       alpaka::wait(queue);
0280 
0281       const GenericSoAConstView& genericViewHostCollection = genericHostCollection.const_view();
0282 
0283       for (size_t i = 0; i < elems; i++) {
0284         REQUIRE(genericViewHostCollection[i].x() == positionViewHostCollection[i].x());
0285         REQUIRE(genericViewHostCollection[i].y() == positionViewHostCollection[i].y());
0286         REQUIRE(genericViewHostCollection[i].z() == positionViewHostCollection[i].z());
0287         REQUIRE(genericViewHostCollection[i].candidateDirection()(0) ==
0288                 pcaViewHostCollection[i].candidateDirection()(0));
0289         REQUIRE(genericViewHostCollection[i].candidateDirection()(1) ==
0290                 pcaViewHostCollection[i].candidateDirection()(1));
0291         REQUIRE(genericViewHostCollection[i].candidateDirection()(2) ==
0292                 pcaViewHostCollection[i].candidateDirection()(2));
0293       }
0294     }
0295   }
0296 }