File indexing completed on 2024-04-06 12:15:41
0001 #include <cassert>
0002 #include <cstdlib>
0003 #include <iostream>
0004 #include <new>
0005
0006 #include <alpaka/alpaka.hpp>
0007
0008 #include "FWCore/Utilities/interface/stringize.h"
0009 #include "HeterogeneousCore/AlpakaInterface/interface/SimpleVector.h"
0010 #include "HeterogeneousCore/AlpakaInterface/interface/devices.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 cms::alpakatools;
0016 using namespace ALPAKA_ACCELERATOR_NAMESPACE;
0017
0018 struct vector_pushback {
0019 template <typename TAcc>
0020 ALPAKA_FN_ACC void operator()(const TAcc& acc, SimpleVector<int>* foo) const {
0021 for (auto index : uniform_elements(acc))
0022 foo->push_back(acc, index);
0023 }
0024 };
0025
0026 struct vector_reset {
0027 template <typename TAcc>
0028 ALPAKA_FN_ACC void operator()(const TAcc& acc, SimpleVector<int>* foo) const {
0029 foo->reset();
0030 }
0031 };
0032
0033 struct vector_emplace_back {
0034 template <typename TAcc>
0035 ALPAKA_FN_ACC void operator()(const TAcc& acc, SimpleVector<int>* foo) const {
0036 for (auto index : uniform_elements(acc))
0037 foo->emplace_back(acc, index);
0038 }
0039 };
0040
0041 int main() {
0042
0043 auto const& devices = cms::alpakatools::devices<Platform>();
0044 if (devices.empty()) {
0045 std::cerr << "No devices available for the " EDM_STRINGIZE(ALPAKA_ACCELERATOR_NAMESPACE) " backend, "
0046 "the test will be skipped.\n";
0047 exit(EXIT_FAILURE);
0048 }
0049
0050
0051 for (auto const& device : devices) {
0052 Queue queue(device);
0053 auto maxN = 10000;
0054 auto vec_h = make_host_buffer<cms::alpakatools::SimpleVector<int>>(queue);
0055 auto vec_d = make_device_buffer<cms::alpakatools::SimpleVector<int>>(queue);
0056 auto data_h = make_host_buffer<int[]>(queue, maxN);
0057 auto data_d = make_device_buffer<int[]>(queue, maxN);
0058
0059 [[maybe_unused]] auto v = make_SimpleVector(maxN, data_d.data());
0060
0061
0062 auto tmp_vec_h = make_host_buffer<cms::alpakatools::SimpleVector<int>>(queue);
0063 make_SimpleVector(tmp_vec_h.data(), maxN, data_d.data());
0064 assert(tmp_vec_h->size() == 0);
0065 assert(tmp_vec_h->capacity() == static_cast<int>(maxN));
0066
0067
0068 alpaka::memcpy(queue, vec_d, tmp_vec_h);
0069 alpaka::wait(queue);
0070
0071 int numBlocks = 5;
0072 int numThreadsPerBlock = 256;
0073 const auto workDiv = make_workdiv<Acc1D>(numBlocks, numThreadsPerBlock);
0074 alpaka::enqueue(queue, alpaka::createTaskKernel<Acc1D>(workDiv, vector_pushback(), vec_d.data()));
0075 alpaka::wait(queue);
0076
0077 alpaka::memcpy(queue, vec_h, vec_d);
0078 alpaka::wait(queue);
0079 printf("vec_h->size()=%d, numBlocks * numThreadsPerBlock=%d, maxN=%d\n",
0080 vec_h->size(),
0081 numBlocks * numThreadsPerBlock,
0082 maxN);
0083 assert(vec_h->size() == (numBlocks * numThreadsPerBlock < maxN ? numBlocks * numThreadsPerBlock : maxN));
0084 alpaka::enqueue(queue, alpaka::createTaskKernel<Acc1D>(workDiv, vector_reset(), vec_d.data()));
0085 alpaka::wait(queue);
0086
0087 alpaka::memcpy(queue, vec_h, vec_d);
0088 alpaka::wait(queue);
0089
0090 assert(vec_h->size() == 0);
0091
0092 alpaka::enqueue(queue, alpaka::createTaskKernel<Acc1D>(workDiv, vector_emplace_back(), vec_d.data()));
0093 alpaka::wait(queue);
0094
0095 alpaka::memcpy(queue, vec_h, vec_d);
0096 alpaka::wait(queue);
0097
0098 assert(vec_h->size() == (numBlocks * numThreadsPerBlock < maxN ? numBlocks * numThreadsPerBlock : maxN));
0099
0100 alpaka::memcpy(queue, data_h, data_d);
0101 }
0102 std::cout << "TEST PASSED" << std::endl;
0103 return 0;
0104 }