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>
70 using size_type = size_t;
81 GpuVector(
const GpuVector<T>& other);
94 explicit GpuVector(
const std::vector<T>& data);
104 GpuVector& operator=(
const GpuVector<T>& other);
117 template<
int BlockDimension>
118 explicit GpuVector(
const Dune::BlockVector<Dune::FieldVector<T, BlockDimension>>& bvector)
119 : GpuVector(bvector.
dim())
121 copyFromHost(bvector);
131 GpuVector& operator=(T scalar);
140 explicit GpuVector(
const size_t numberOfElements);
154 GpuVector(
const T* dataOnHost,
const size_t numberOfElements);
159 virtual ~GpuVector();
169 const T* data()
const;
178 template <
int BlockDimension>
179 void copyFromHost(
const Dune::BlockVector<Dune::FieldVector<T, BlockDimension>>& bvector)
183 OPM_THROW(std::runtime_error,
184 fmt::format(
"Given incompatible vector size. GpuVector has size {}, \n"
185 "however, BlockVector has N() = {}, and dim = {}.",
190 const auto dataPointer =
static_cast<const T*
>(&(bvector[0][0]));
191 copyFromHost(dataPointer, m_numberOfElements);
202 template <
int BlockDimension>
203 void copyFromHostAsync(
const Dune::BlockVector<Dune::FieldVector<T, BlockDimension>>& bvector, cudaStream_t stream =
detail::DEFAULT_STREAM)
207 OPM_THROW(std::runtime_error,
208 fmt::format(
"Given incompatible vector size. GpuVector has size {}, \n"
209 "however, BlockVector has N() = {}, and dim = {}.",
214 const auto dataPointer =
static_cast<const T*
>(&(bvector[0][0]));
215 copyFromHostAsync(dataPointer, m_numberOfElements, stream);
225 template <
int BlockDimension>
226 void copyToHost(Dune::BlockVector<Dune::FieldVector<T, BlockDimension>>& bvector)
const
230 OPM_THROW(std::runtime_error,
231 fmt::format(
"Given incompatible vector size. GpuVector has size {},\n however, the BlockVector "
232 "has has N() = {}, and dim() = {}.",
237 const auto dataPointer =
static_cast<T*
>(&(bvector[0][0]));
238 copyToHost(dataPointer, m_numberOfElements);
249 template <
int BlockDimension>
250 void copyToHostAsync(Dune::BlockVector<Dune::FieldVector<T, BlockDimension>>& bvector, cudaStream_t stream =
detail::DEFAULT_STREAM)
const
254 OPM_THROW(std::runtime_error,
255 fmt::format(
"Given incompatible vector size. GpuVector has size {},\n however, the BlockVector "
256 "has has N() = {}, and dim() = {}.",
261 const auto dataPointer =
static_cast<T*
>(&(bvector[0][0]));
262 copyToHostAsync(dataPointer, m_numberOfElements, stream);
272 void copyFromHost(
const T* dataPointer,
size_t numberOfElements);
284 void copyFromHostAsync(
const T* dataPointer,
size_t numberOfElements, cudaStream_t stream =
detail::DEFAULT_STREAM);
295 void copyToHost(T* dataPointer,
size_t numberOfElements)
const;
307 void copyToHostAsync(T* dataPointer,
size_t numberOfElements, cudaStream_t stream =
detail::DEFAULT_STREAM)
const;
318 void copyFromHost(
const std::vector<T>& data);
337 void copyToHost(std::vector<T>& data)
const;
356 void copyFromDeviceToDevice(
const GpuVector<T>& other)
const;
358 void prepareSendBuf(GpuVector<T>& buffer,
const GpuVector<int>& indexSet)
const;
359 void syncFromRecvBuf(GpuVector<T>& buffer,
const GpuVector<int>& indexSet)
const;
369 GpuVector<T>& operator*=(
const T& scalar);
379 GpuVector<T>& axpy(T alpha,
const GpuVector<T>& y);
387 GpuVector<T>& operator+=(
const GpuVector<T>& other);
395 GpuVector<T>& operator-=(
const GpuVector<T>& other);
405 T dot(
const GpuVector<T>& other)
const;
421 T dot(
const GpuVector<T>& other,
const GpuVector<int>& indexSet, GpuVector<T>& buffer)
const;
428 T two_norm(
const GpuVector<int>& indexSet, GpuVector<T>& buffer)
const;
436 T dot(
const GpuVector<T>& other,
const GpuVector<int>& indexSet)
const;
443 T two_norm(
const GpuVector<int>& indexSet)
const;
450 size_type
dim()
const;
457 std::vector<T> asStdVector()
const;
463 template <
int blockSize>
464 Dune::BlockVector<Dune::FieldVector<T, blockSize>> asDuneBlockVector()
const
466 OPM_ERROR_IF(
dim() % blockSize != 0,
467 fmt::format(
"blockSize is not a multiple of dim(). Given blockSize = {}, and dim() = {}",
471 Dune::BlockVector<Dune::FieldVector<T, blockSize>> returnValue(
dim() / blockSize);
472 copyToHost(returnValue);
495 std::vector<T> v = asStdVector();
496 std::string res =
"";
505 T* m_dataOnDevice =
nullptr;
509 const int m_numberOfElements;
510 detail::CuBlasHandle& m_cuBlasHandle;
512 void assertSameSize(
const GpuVector<T>& other)
const;
513 void assertSameSize(
int size)
const;
515 void assertHasElements()
const;
static constexpr int dim
Definition: structuredgridvanguard.hh:68
void syncFromRecvBuf(T *deviceA, T *buffer, size_t numberOfElements, const int *indices)
constexpr cudaStream_t DEFAULT_STREAM
The default GPU stream (stream 0)
Definition: gpu_constants.hpp:31
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: AmgxInterface.hpp:38
std::string toDebugString()
Definition: GpuVector.hpp:493
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)