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));
98 GpuBuffer(GpuBuffer<T>&& other) noexcept
99 : m_dataOnDevice(other.m_dataOnDevice)
100 , m_numberOfElements(other.m_numberOfElements)
102 other.m_dataOnDevice =
nullptr;
103 other.m_numberOfElements = 0;
114 GpuBuffer<T>& operator=(GpuBuffer<T>&& other)
noexcept
116 if (
this != &other) {
118 m_dataOnDevice = other.m_dataOnDevice;
119 m_numberOfElements = other.m_numberOfElements;
120 other.m_dataOnDevice =
nullptr;
121 other.m_numberOfElements = 0;
135 explicit GpuBuffer(
const std::vector<T>& data)
136 : GpuBuffer(data.size())
146 explicit GpuBuffer(
const size_t numberOfElements)
147 : m_numberOfElements(numberOfElements)
162 GpuBuffer(
const T* dataOnHost,
const size_t numberOfElements)
163 : GpuBuffer(numberOfElements)
166 m_dataOnDevice, dataOnHost, m_numberOfElements *
sizeof(T), cudaMemcpyHostToDevice));
183 return m_dataOnDevice;
189 const T* data()
const
191 return m_dataOnDevice;
201 template <
int BlockDimension>
202 void copyFromHost(
const Dune::BlockVector<Dune::FieldVector<T, BlockDimension>>& bvector)
205 if (m_numberOfElements != bvector.size()) {
206 OPM_THROW(std::runtime_error,
207 fmt::format(
"Given incompatible vector size. GpuBuffer has size {}, \n"
208 "however, BlockVector has N() = {}, and size = {}.",
213 const auto dataPointer =
static_cast<const T*
>(&(bvector[0][0]));
214 copyFromHost(dataPointer, m_numberOfElements);
224 template <
int BlockDimension>
225 void copyToHost(Dune::BlockVector<Dune::FieldVector<T, BlockDimension>>& bvector)
const
228 if (m_numberOfElements != bvector.size()) {
229 OPM_THROW(std::runtime_error,
230 fmt::format(
"Given incompatible vector size. GpuBuffer has size {},\n however, the BlockVector "
231 "has has N() = {}, and size() = {}.",
236 const auto dataPointer =
static_cast<T*
>(&(bvector[0][0]));
237 copyToHost(dataPointer, m_numberOfElements);
247 void copyFromHost(
const T* dataPointer,
size_t numberOfElements)
249 if (numberOfElements > size()) {
250 OPM_THROW(std::runtime_error,
251 fmt::format(fmt::runtime(
"Requesting to copy too many elements. "
252 "buffer has {} elements, while {} was requested."),
256 OPM_GPU_SAFE_CALL(cudaMemcpy(data(), dataPointer, numberOfElements *
sizeof(T), cudaMemcpyHostToDevice));
266 void copyToHost(T* dataPointer,
size_t numberOfElements)
const
268 assertSameSize(numberOfElements);
269 OPM_GPU_SAFE_CALL(cudaMemcpy(dataPointer, data(), numberOfElements *
sizeof(T), cudaMemcpyDeviceToHost));
279 void copyFromHost(
const std::vector<T>& data)
281 assertSameSize(data.size());
287 if constexpr (std::is_same_v<T, bool>)
289 auto tmp = std::make_unique<bool[]>(data.size());
290 for (
size_t i = 0; i < data.size(); ++i) {
291 tmp[i] =
static_cast<bool>(data[i]);
293 copyFromHost(tmp.get(), data.size());
296 copyFromHost(data.data(), data.size());
307 void copyToHost(std::vector<T>& data)
const
309 assertSameSize(data.size());
315 if constexpr (std::is_same_v<T, bool>)
317 auto tmp = std::make_unique<bool[]>(data.size());
318 copyToHost(tmp.get(), data.size());
319 for (
size_t i = 0; i < data.size(); ++i) {
320 data[i] =
static_cast<bool>(tmp[i]);
325 copyToHost(data.data(), data.size());
333 size_type size()
const
335 return m_numberOfElements;
342 void resize(
size_t newSize)
345 OPM_THROW(std::invalid_argument,
"Setting a GpuBuffer size to a non-positive number is not allowed");
348 if (m_numberOfElements == 0) {
354 T* tmpBuffer =
nullptr;
358 size_t sizeOfMove = std::min({m_numberOfElements, newSize});
361 sizeOfMove *
sizeof(T),
362 cudaMemcpyDeviceToDevice));
368 m_dataOnDevice = tmpBuffer;
372 m_numberOfElements = newSize;
379 std::vector<T> asStdVector()
const
381 std::vector<T> temporary(m_numberOfElements);
382 copyToHost(temporary);
387 T* m_dataOnDevice =
nullptr;
388 size_t m_numberOfElements = 0;
390 void assertSameSize(
const GpuBuffer<T>& other)
const
392 assertSameSize(other.m_numberOfElements);
395 void assertSameSize(
size_t size)
const
397 if (size != m_numberOfElements) {
398 OPM_THROW(std::invalid_argument,
399 fmt::format(fmt::runtime(
"Given buffer has {}, while we have {}."),
400 size, m_numberOfElements));
404 void assertHasElements()
const
406 if (m_numberOfElements <= 0) {
407 OPM_THROW(std::invalid_argument,
"We have 0 elements");
414 return GpuView<T>(buf.data(), buf.size());
418GpuView<const T>
make_view(
const GpuBuffer<T>& buf) {
419 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:164
#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:185
A small, fixed‑dimension MiniVector class backed by std::array that can be used in both host and CUDA...
Definition: ThermalGasWaterFlowProblem.hpp:95
ThermalGasWaterFlowProblem< Scalar, GpuView > make_view(ThermalGasWaterFlowProblem< Scalar, GpuBuffer > &buffer)
Definition: ThermalGasWaterFlowProblem.hpp:111