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>
26#include <opm/simulators/linalg/gpuistl/detail/CuBlasHandle.hpp>
27#include <opm/simulators/linalg/gpuistl/detail/safe_conversion.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);
112 GpuVector& operator=(T scalar);
121 explicit GpuVector(
const size_t numberOfElements);
135 GpuVector(
const T* dataOnHost,
const size_t numberOfElements);
140 virtual ~GpuVector();
150 const T* data()
const;
159 template <
int BlockDimension>
160 void copyFromHost(
const Dune::BlockVector<Dune::FieldVector<T, BlockDimension>>& bvector)
163 if (detail::to_size_t(m_numberOfElements) != bvector.dim()) {
164 OPM_THROW(std::runtime_error,
165 fmt::format(
"Given incompatible vector size. GpuVector has size {}, \n"
166 "however, BlockVector has N() = {}, and dim = {}.",
171 const auto dataPointer =
static_cast<const T*
>(&(bvector[0][0]));
172 copyFromHost(dataPointer, m_numberOfElements);
182 template <
int BlockDimension>
183 void copyToHost(Dune::BlockVector<Dune::FieldVector<T, BlockDimension>>& bvector)
const
186 if (detail::to_size_t(m_numberOfElements) != bvector.dim()) {
187 OPM_THROW(std::runtime_error,
188 fmt::format(
"Given incompatible vector size. GpuVector has size {},\n however, the BlockVector "
189 "has has N() = {}, and dim() = {}.",
194 const auto dataPointer =
static_cast<T*
>(&(bvector[0][0]));
195 copyToHost(dataPointer, m_numberOfElements);
205 void copyFromHost(
const T* dataPointer,
size_t numberOfElements);
206 void copyFromHost(
const T* dataPointer,
size_t numberOfElements, cudaStream_t stream);
215 void copyToHost(T* dataPointer,
size_t numberOfElements)
const;
224 void copyFromHost(
const std::vector<T>& data);
233 void copyToHost(std::vector<T>& data)
const;
235 void prepareSendBuf(GpuVector<T>& buffer,
const GpuVector<int>& indexSet)
const;
236 void syncFromRecvBuf(GpuVector<T>& buffer,
const GpuVector<int>& indexSet)
const;
246 GpuVector<T>& operator*=(
const T& scalar);
256 GpuVector<T>& axpy(T alpha,
const GpuVector<T>& y);
264 GpuVector<T>& operator+=(
const GpuVector<T>& other);
272 GpuVector<T>& operator-=(
const GpuVector<T>& other);
282 T dot(
const GpuVector<T>& other)
const;
298 T dot(
const GpuVector<T>& other,
const GpuVector<int>& indexSet, GpuVector<T>& buffer)
const;
305 T two_norm(
const GpuVector<int>& indexSet, GpuVector<T>& buffer)
const;
313 T dot(
const GpuVector<T>& other,
const GpuVector<int>& indexSet)
const;
320 T two_norm(
const GpuVector<int>& indexSet)
const;
327 size_type dim()
const;
334 std::vector<T> asStdVector()
const;
340 template <
int blockSize>
341 Dune::BlockVector<Dune::FieldVector<T, blockSize>> asDuneBlockVector()
const
343 OPM_ERROR_IF(dim() % blockSize != 0,
344 fmt::format(
"blockSize is not a multiple of dim(). Given blockSize = {}, and dim() = {}",
348 Dune::BlockVector<Dune::FieldVector<T, blockSize>> returnValue(dim() / blockSize);
349 copyToHost(returnValue);
370 std::string toDebugString()
372 std::vector<T> v = asStdVector();
373 std::string res =
"";
375 res += std::to_string(element) +
" ";
377 res += std::to_string(v[v.size()-1]);
382 T* m_dataOnDevice =
nullptr;
386 const int m_numberOfElements;
387 detail::CuBlasHandle& m_cuBlasHandle;
389 void assertSameSize(
const GpuVector<T>& other)
const;
390 void assertSameSize(
int size)
const;
392 void assertHasElements()
const;
void setZeroAtIndexSet(T *deviceData, size_t numberOfElements, const int *indices)
setZeroAtIndexSet sets deviceData to zero in the indices of contained in indices
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:242