19#ifndef OPM_CUVECTOR_HEADER_HPP
20#define OPM_CUVECTOR_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/cuistl/detail/CuBlasHandle.hpp>
27#include <opm/simulators/linalg/cuistl/detail/safe_conversion.hpp>
68 using size_type = size_t;
79 CuVector(
const CuVector<T>& other);
92 explicit CuVector(
const std::vector<T>& data);
102 CuVector& operator=(
const CuVector<T>& other);
111 CuVector& operator=(T scalar);
120 explicit CuVector(
const size_t numberOfElements);
134 CuVector(
const T* dataOnHost,
const size_t numberOfElements);
149 const T* data()
const;
158 template <
int BlockDimension>
159 void copyFromHost(
const Dune::BlockVector<Dune::FieldVector<T, BlockDimension>>& bvector)
162 if (detail::to_size_t(m_numberOfElements) != bvector.dim()) {
163 OPM_THROW(std::runtime_error,
164 fmt::format(
"Given incompatible vector size. CuVector has size {}, \n"
165 "however, BlockVector has N() = {}, and dim = {}.",
170 const auto dataPointer =
static_cast<const T*
>(&(bvector[0][0]));
171 copyFromHost(dataPointer, m_numberOfElements);
181 template <
int BlockDimension>
182 void copyToHost(Dune::BlockVector<Dune::FieldVector<T, BlockDimension>>& bvector)
const
185 if (detail::to_size_t(m_numberOfElements) != bvector.dim()) {
186 OPM_THROW(std::runtime_error,
187 fmt::format(
"Given incompatible vector size. CuVector has size {},\n however, the BlockVector "
188 "has has N() = {}, and dim() = {}.",
193 const auto dataPointer =
static_cast<T*
>(&(bvector[0][0]));
194 copyToHost(dataPointer, m_numberOfElements);
204 void copyFromHost(
const T* dataPointer,
size_t numberOfElements);
213 void copyToHost(T* dataPointer,
size_t numberOfElements)
const;
222 void copyFromHost(
const std::vector<T>& data);
231 void copyToHost(std::vector<T>& data)
const;
241 CuVector<T>& operator*=(
const T& scalar);
251 CuVector<T>& axpy(T alpha,
const CuVector<T>& y);
259 CuVector<T>& operator+=(
const CuVector<T>& other);
267 CuVector<T>& operator-=(
const CuVector<T>& other);
277 T dot(
const CuVector<T>& other)
const;
293 T dot(
const CuVector<T>& other,
const CuVector<int>& indexSet, CuVector<T>& buffer)
const;
300 T two_norm(
const CuVector<int>& indexSet, CuVector<T>& buffer)
const;
308 T dot(
const CuVector<T>& other,
const CuVector<int>& indexSet)
const;
315 T two_norm(
const CuVector<int>& indexSet)
const;
322 size_type dim()
const;
329 std::vector<T> asStdVector()
const;
335 template <
int blockSize>
336 Dune::BlockVector<Dune::FieldVector<T, blockSize>> asDuneBlockVector()
const
338 OPM_ERROR_IF(dim() % blockSize != 0,
339 fmt::format(
"blockSize is not a multiple of dim(). Given blockSize = {}, and dim() = {}",
343 Dune::BlockVector<Dune::FieldVector<T, blockSize>> returnValue(dim() / blockSize);
344 copyToHost(returnValue);
362 void setZeroAtIndexSet(
const CuVector<int>& indexSet);
365 T* m_dataOnDevice =
nullptr;
369 const int m_numberOfElements;
370 detail::CuBlasHandle& m_cuBlasHandle;
372 void assertSameSize(
const CuVector<T>& other)
const;
373 void assertSameSize(
int size)
const;
375 void assertHasElements()
const;