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>
59 using size_type = size_t;
70 GpuBuffer(
const GpuBuffer<T>& other);
81 explicit GpuBuffer(
const std::vector<T>& data);
93 explicit GpuBuffer(
const size_t numberOfElements);
105 GpuBuffer(
const T* dataOnHost,
const size_t numberOfElements);
110 virtual ~GpuBuffer();
120 const T* data()
const;
125 __host__ __device__ T& front()
130 return m_dataOnDevice[0];
136 __host__ __device__ T& back()
141 return m_dataOnDevice[m_numberOfElements-1];
147 __host__ __device__ T front()
const
152 return m_dataOnDevice[0];
158 __host__ __device__ T back()
const
163 return m_dataOnDevice[m_numberOfElements-1];
173 template <
int BlockDimension>
174 void copyFromHost(
const Dune::BlockVector<Dune::FieldVector<T, BlockDimension>>& bvector)
177 if (m_numberOfElements != bvector.size()) {
178 OPM_THROW(std::runtime_error,
179 fmt::format(
"Given incompatible vector size. GpuBuffer has size {}, \n"
180 "however, BlockVector has N() = {}, and size = {}.",
185 const auto dataPointer =
static_cast<const T*
>(&(bvector[0][0]));
186 copyFromHost(dataPointer, m_numberOfElements);
196 template <
int BlockDimension>
197 void copyToHost(Dune::BlockVector<Dune::FieldVector<T, BlockDimension>>& bvector)
const
200 if (m_numberOfElements != bvector.size()) {
201 OPM_THROW(std::runtime_error,
202 fmt::format(
"Given incompatible vector size. GpuBuffer has size {},\n however, the BlockVector "
203 "has has N() = {}, and size() = {}.",
208 const auto dataPointer =
static_cast<T*
>(&(bvector[0][0]));
209 copyToHost(dataPointer, m_numberOfElements);
219 void copyFromHost(
const T* dataPointer,
size_t numberOfElements);
228 void copyToHost(T* dataPointer,
size_t numberOfElements)
const;
237 void copyFromHost(
const std::vector<T>& data);
246 void copyToHost(std::vector<T>& data)
const;
252 size_type size()
const;
264 std::vector<T> asStdVector()
const;
267 T* m_dataOnDevice =
nullptr;
268 size_t m_numberOfElements;
270 void assertSameSize(
const GpuBuffer<T>& other)
const;
271 void assertSameSize(
size_t size)
const;
273 void assertHasElements()
const;
277GpuView<const T> make_view(
const GpuBuffer<T>&);
Definition: autotuner.hpp:29