File indexing completed on 2024-11-11 23:31:37
0001 #ifndef HeterogeneousCore_AlpakaTest_interface_AlpakaESTestData_h
0002 #define HeterogeneousCore_AlpakaTest_interface_AlpakaESTestData_h
0003
0004 #include "DataFormats/Portable/interface/PortableHostCollection.h"
0005 #include "DataFormats/Portable/interface/PortableCollection.h"
0006 #include "HeterogeneousCore/AlpakaInterface/interface/CopyToDevice.h"
0007 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0008 #include "HeterogeneousCore/AlpakaInterface/interface/memory.h"
0009 #include "HeterogeneousCore/AlpakaTest/interface/AlpakaESTestSoA.h"
0010
0011 namespace cms::alpakatest {
0012
0013 using AlpakaESTestDataAHost = PortableHostCollection<AlpakaESTestSoAA>;
0014 using AlpakaESTestDataCHost = PortableHostCollection<AlpakaESTestSoAC>;
0015 using AlpakaESTestDataDHost = PortableHostCollection<AlpakaESTestSoAD>;
0016
0017 using AlpakaESTestDataACMultiHost = PortableHostMultiCollection<AlpakaESTestSoAA, AlpakaESTestSoAC>;
0018
0019
0020 template <typename TDev>
0021 class AlpakaESTestDataB {
0022 public:
0023 using Buffer = cms::alpakatools::device_buffer<TDev, int[]>;
0024 using ConstBuffer = cms::alpakatools::const_device_buffer<TDev, int[]>;
0025
0026 explicit AlpakaESTestDataB(Buffer buffer) : buffer_(std::move(buffer)) {}
0027
0028 Buffer buffer() { return buffer_; }
0029 ConstBuffer buffer() const { return buffer_; }
0030 ConstBuffer const_buffer() const { return buffer_; }
0031
0032 int const* data() const { return buffer_.data(); }
0033 auto size() const { return alpaka::getExtentProduct(buffer_); }
0034
0035 private:
0036 Buffer buffer_;
0037 };
0038
0039
0040
0041 template <typename TDev>
0042 class AlpakaESTestDataE {
0043 public:
0044 using ECollection = PortableCollection<AlpakaESTestSoAE, TDev>;
0045 using EDataCollection = PortableCollection<AlpakaESTestSoAEData, TDev>;
0046
0047 class ConstView {
0048 public:
0049 constexpr ConstView(typename ECollection::ConstView e, typename EDataCollection::ConstView data)
0050 : eView_(e), dataView_(data) {}
0051
0052 constexpr auto size() const { return eView_.metadata().size(); }
0053 constexpr int val(int i) const { return eView_.val(i); }
0054 constexpr int val2(int i) const { return dataView_.val2(eView_.ind(i)); }
0055
0056 private:
0057 typename ECollection::ConstView eView_;
0058 typename EDataCollection::ConstView dataView_;
0059 };
0060
0061 AlpakaESTestDataE(size_t size, size_t dataSize) : e_(size), data_(dataSize) {}
0062
0063 AlpakaESTestDataE(ECollection e, EDataCollection data) : e_(std::move(e)), data_(std::move(data)) {}
0064
0065 ECollection const& e() const { return e_; }
0066 EDataCollection const& data() const { return data_; }
0067
0068 ConstView view() const { return const_view(); }
0069 ConstView const_view() const { return ConstView(e_.const_view(), data_.const_view()); }
0070
0071 private:
0072 ECollection e_;
0073 EDataCollection data_;
0074 };
0075 using AlpakaESTestDataEHost = AlpakaESTestDataE<alpaka_common::DevHost>;
0076
0077 }
0078
0079 namespace cms::alpakatools {
0080
0081
0082
0083
0084 template <>
0085 struct CopyToDevice<cms::alpakatest::AlpakaESTestDataB<alpaka_common::DevHost>> {
0086 template <typename TQueue>
0087 static auto copyAsync(TQueue& queue, cms::alpakatest::AlpakaESTestDataB<alpaka_common::DevHost> const& srcData) {
0088
0089
0090
0091
0092
0093
0094
0095
0096
0097
0098 auto dstBuffer = cms::alpakatools::make_device_buffer<int[]>(queue, srcData.size());
0099 alpaka::memcpy(queue, dstBuffer, srcData.buffer());
0100 return cms::alpakatest::AlpakaESTestDataB<alpaka::Dev<TQueue>>(std::move(dstBuffer));
0101 }
0102 };
0103
0104 template <>
0105 struct CopyToDevice<cms::alpakatest::AlpakaESTestDataEHost> {
0106 template <typename TQueue>
0107 static auto copyAsync(TQueue& queue, cms::alpakatest::AlpakaESTestDataEHost const& srcData) {
0108 using ECopy = CopyToDevice<cms::alpakatest::AlpakaESTestDataEHost::ECollection>;
0109 using EDataCopy = CopyToDevice<cms::alpakatest::AlpakaESTestDataEHost::EDataCollection>;
0110 using TDevice = alpaka::Dev<TQueue>;
0111 return cms::alpakatest::AlpakaESTestDataE<TDevice>(ECopy::copyAsync(queue, srcData.e()),
0112 EDataCopy::copyAsync(queue, srcData.data()));
0113 }
0114 };
0115 }
0116
0117 #endif