19#ifndef OPM_GPUBUFFER_HEADER_HPP
20#define OPM_GPUBUFFER_HEADER_HPP
21#include <dune/common/fvector.hh>
22#include <dune/istl/bvector.hh>
25#include <opm/common/ErrorMacros.hpp>
31#include <cuda_runtime.h>
61 using size_type = size_t;
67 GpuBuffer() =
default;
77 GpuBuffer(
const GpuBuffer<T>& other)
78 : GpuBuffer(other.m_numberOfElements)
80 assertSameSize(other);
81 if (m_numberOfElements == 0) {
86 m_numberOfElements *
sizeof(T),
87 cudaMemcpyDeviceToDevice));
99 explicit GpuBuffer(
const std::vector<T>& data)
100 : GpuBuffer(data.size())
110 explicit GpuBuffer(
const size_t numberOfElements)
111 : m_numberOfElements(numberOfElements)
126 GpuBuffer(
const T* dataOnHost,
const size_t numberOfElements)
127 : GpuBuffer(numberOfElements)
130 m_dataOnDevice, dataOnHost, m_numberOfElements *
sizeof(T), cudaMemcpyHostToDevice));
147 return m_dataOnDevice;
153 const T* data()
const
155 return m_dataOnDevice;
165 template <
int BlockDimension>
166 void copyFromHost(
const Dune::BlockVector<Dune::FieldVector<T, BlockDimension>>& bvector)
169 if (m_numberOfElements != bvector.size()) {
170 OPM_THROW(std::runtime_error,
171 fmt::format(
"Given incompatible vector size. GpuBuffer has size {}, \n"
172 "however, BlockVector has N() = {}, and size = {}.",
177 const auto dataPointer =
static_cast<const T*
>(&(bvector[0][0]));
178 copyFromHost(dataPointer, m_numberOfElements);
188 template <
int BlockDimension>
189 void copyToHost(Dune::BlockVector<Dune::FieldVector<T, BlockDimension>>& bvector)
const
192 if (m_numberOfElements != bvector.size()) {
193 OPM_THROW(std::runtime_error,
194 fmt::format(
"Given incompatible vector size. GpuBuffer has size {},\n however, the BlockVector "
195 "has has N() = {}, and size() = {}.",
200 const auto dataPointer =
static_cast<T*
>(&(bvector[0][0]));
201 copyToHost(dataPointer, m_numberOfElements);
211 void copyFromHost(
const T* dataPointer,
size_t numberOfElements)
213 if (numberOfElements > size()) {
214 OPM_THROW(std::runtime_error,
215 fmt::format(
"Requesting to copy too many elements. buffer has {} elements, while {} was requested.",
219 OPM_GPU_SAFE_CALL(cudaMemcpy(data(), dataPointer, numberOfElements *
sizeof(T), cudaMemcpyHostToDevice));
229 void copyToHost(T* dataPointer,
size_t numberOfElements)
const
231 assertSameSize(numberOfElements);
232 OPM_GPU_SAFE_CALL(cudaMemcpy(dataPointer, data(), numberOfElements *
sizeof(T), cudaMemcpyDeviceToHost));
242 void copyFromHost(
const std::vector<T>& data)
244 assertSameSize(data.size());
250 if constexpr (std::is_same_v<T, bool>)
252 auto tmp = std::make_unique<bool[]>(data.size());
253 for (
size_t i = 0; i < data.size(); ++i) {
254 tmp[i] =
static_cast<bool>(data[i]);
256 copyFromHost(tmp.get(), data.size());
259 copyFromHost(data.data(), data.size());
270 void copyToHost(std::vector<T>& data)
const
272 assertSameSize(data.size());
278 if constexpr (std::is_same_v<T, bool>)
280 auto tmp = std::make_unique<bool[]>(data.size());
281 copyToHost(tmp.get(), data.size());
282 for (
size_t i = 0; i < data.size(); ++i) {
283 data[i] =
static_cast<bool>(tmp[i]);
288 copyToHost(data.data(), data.size());
296 size_type size()
const
298 return m_numberOfElements;
305 void resize(
size_t newSize)
308 OPM_THROW(std::invalid_argument,
"Setting a GpuBuffer size to a non-positive number is not allowed");
311 if (m_numberOfElements == 0) {
317 T* tmpBuffer =
nullptr;
321 size_t sizeOfMove = std::min({m_numberOfElements, newSize});
324 sizeOfMove *
sizeof(T),
325 cudaMemcpyDeviceToDevice));
331 m_dataOnDevice = tmpBuffer;
335 m_numberOfElements = newSize;
342 std::vector<T> asStdVector()
const
344 std::vector<T> temporary(m_numberOfElements);
345 copyToHost(temporary);
350 T* m_dataOnDevice =
nullptr;
351 size_t m_numberOfElements = 0;
353 void assertSameSize(
const GpuBuffer<T>& other)
const
355 assertSameSize(other.m_numberOfElements);
358 void assertSameSize(
size_t size)
const
360 if (size != m_numberOfElements) {
361 OPM_THROW(std::invalid_argument,
362 fmt::format(
"Given buffer has {}, while we have {}.", size, m_numberOfElements));
366 void assertHasElements()
const
368 if (m_numberOfElements <= 0) {
369 OPM_THROW(std::invalid_argument,
"We have 0 elements");
376 return GpuView<T>(buf.data(), buf.size());
380GpuView<const T>
make_view(
const GpuBuffer<T>& buf) {
381 return GpuView<const T>(buf.data(), buf.size());
#define OPM_GPU_SAFE_CALL(expression)
OPM_GPU_SAFE_CALL checks the return type of the GPU expression (function call) and throws an exceptio...
Definition: gpu_safe_call.hpp:150
#define OPM_GPU_WARN_IF_ERROR(expression)
OPM_GPU_WARN_IF_ERROR checks the return type of the GPU expression (function call) and issues a warni...
Definition: gpu_safe_call.hpp:171
A small, fixed‑dimension MiniVector class backed by std::array that can be used in both host and CUDA...
Definition: AmgxInterface.hpp:38
PointerView< T > make_view(const std::shared_ptr< T > &ptr)
Definition: gpu_smart_pointer.hpp:318