19#ifndef OPM_GPUVECTOR_HEADER_HPP
20#define OPM_GPUVECTOR_HEADER_HPP
21#include <dune/common/fvector.hh>
22#include <dune/istl/bvector.hh>
25#include <opm/common/ErrorMacros.hpp>
69 using size_type = size_t;
80 GpuVector(
const GpuVector<T>& other);
93 explicit GpuVector(
const std::vector<T>& data);
103 GpuVector& operator=(
const GpuVector<T>& other);
116 template<
int BlockDimension>
117 explicit GpuVector(
const Dune::BlockVector<Dune::FieldVector<T, BlockDimension>>& bvector)
118 : GpuVector(bvector.
dim())
120 copyFromHost(bvector);
130 GpuVector& operator=(T scalar);
139 explicit GpuVector(
const size_t numberOfElements);
153 GpuVector(
const T* dataOnHost,
const size_t numberOfElements);
158 virtual ~GpuVector();
168 const T* data()
const;
177 template <
int BlockDimension>
178 void copyFromHost(
const Dune::BlockVector<Dune::FieldVector<T, BlockDimension>>& bvector)
182 OPM_THROW(std::runtime_error,
183 fmt::format(
"Given incompatible vector size. GpuVector has size {}, \n"
184 "however, BlockVector has N() = {}, and dim = {}.",
189 const auto dataPointer =
static_cast<const T*
>(&(bvector[0][0]));
190 copyFromHost(dataPointer, m_numberOfElements);
200 template <
int BlockDimension>
201 void copyToHost(Dune::BlockVector<Dune::FieldVector<T, BlockDimension>>& bvector)
const
205 OPM_THROW(std::runtime_error,
206 fmt::format(
"Given incompatible vector size. GpuVector has size {},\n however, the BlockVector "
207 "has has N() = {}, and dim() = {}.",
212 const auto dataPointer =
static_cast<T*
>(&(bvector[0][0]));
213 copyToHost(dataPointer, m_numberOfElements);
223 void copyFromHost(
const T* dataPointer,
size_t numberOfElements);
224 void copyFromHost(
const T* dataPointer,
size_t numberOfElements, cudaStream_t stream);
233 void copyToHost(T* dataPointer,
size_t numberOfElements)
const;
242 void copyFromHost(
const std::vector<T>& data);
251 void copyToHost(std::vector<T>& data)
const;
257 void copyFromDeviceToDevice(
const GpuVector<T>& other)
const;
259 void prepareSendBuf(GpuVector<T>& buffer,
const GpuVector<int>& indexSet)
const;
260 void syncFromRecvBuf(GpuVector<T>& buffer,
const GpuVector<int>& indexSet)
const;
270 GpuVector<T>& operator*=(
const T& scalar);
280 GpuVector<T>& axpy(T alpha,
const GpuVector<T>& y);
288 GpuVector<T>& operator+=(
const GpuVector<T>& other);
296 GpuVector<T>& operator-=(
const GpuVector<T>& other);
306 T dot(
const GpuVector<T>& other)
const;
322 T dot(
const GpuVector<T>& other,
const GpuVector<int>& indexSet, GpuVector<T>& buffer)
const;
329 T two_norm(
const GpuVector<int>& indexSet, GpuVector<T>& buffer)
const;
337 T dot(
const GpuVector<T>& other,
const GpuVector<int>& indexSet)
const;
344 T two_norm(
const GpuVector<int>& indexSet)
const;
351 size_type
dim()
const;
358 std::vector<T> asStdVector()
const;
364 template <
int blockSize>
365 Dune::BlockVector<Dune::FieldVector<T, blockSize>> asDuneBlockVector()
const
367 OPM_ERROR_IF(
dim() % blockSize != 0,
368 fmt::format(
"blockSize is not a multiple of dim(). Given blockSize = {}, and dim() = {}",
372 Dune::BlockVector<Dune::FieldVector<T, blockSize>> returnValue(
dim() / blockSize);
373 copyToHost(returnValue);
396 std::vector<T> v = asStdVector();
397 std::string res =
"";
406 T* m_dataOnDevice =
nullptr;
410 const int m_numberOfElements;
411 detail::CuBlasHandle& m_cuBlasHandle;
413 void assertSameSize(
const GpuVector<T>& other)
const;
414 void assertSameSize(
int size)
const;
416 void assertHasElements()
const;
static constexpr int dim
Definition: structuredgridvanguard.hh:68
void syncFromRecvBuf(T *deviceA, T *buffer, size_t numberOfElements, const int *indices)
void prepareSendBuf(const T *deviceA, T *buffer, size_t numberOfElements, const int *indices)
__host__ __device__ std::size_t to_size_t(int i)
to_size_t converts a (on most relevant platforms) a 32 bit signed int to a 64 bits unsigned int
Definition: safe_conversion.hpp:86
Definition: autotuner.hpp:30
std::string toDebugString()
Definition: GpuVector.hpp:394
void setZeroAtIndexSet(const GpuVector< int > &indexSet)
The GpuVector class is a simple (arithmetic) vector class for the GPU.
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)