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>
69 using size_type = size_t;
80 CuVector(
const CuVector<T>& other);
93 explicit CuVector(
const std::vector<T>& data);
103 CuVector& operator=(
const CuVector<T>& other);
112 CuVector& operator=(T scalar);
121 explicit CuVector(
const size_t numberOfElements);
135 CuVector(
const T* dataOnHost,
const size_t numberOfElements);
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. CuVector 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. CuVector 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);
214 void copyToHost(T* dataPointer,
size_t numberOfElements)
const;
223 void copyFromHost(
const std::vector<T>& data);
232 void copyToHost(std::vector<T>& data)
const;
242 CuVector<T>& operator*=(
const T& scalar);
252 CuVector<T>& axpy(T alpha,
const CuVector<T>& y);
260 CuVector<T>& operator+=(
const CuVector<T>& other);
268 CuVector<T>& operator-=(
const CuVector<T>& other);
278 T dot(
const CuVector<T>& other)
const;
294 T dot(
const CuVector<T>& other,
const CuVector<int>& indexSet, CuVector<T>& buffer)
const;
301 T two_norm(
const CuVector<int>& indexSet, CuVector<T>& buffer)
const;
309 T dot(
const CuVector<T>& other,
const CuVector<int>& indexSet)
const;
316 T two_norm(
const CuVector<int>& indexSet)
const;
323 size_type dim()
const;
330 std::vector<T> asStdVector()
const;
336 template <
int blockSize>
337 Dune::BlockVector<Dune::FieldVector<T, blockSize>> asDuneBlockVector()
const
339 OPM_ERROR_IF(dim() % blockSize != 0,
340 fmt::format(
"blockSize is not a multiple of dim(). Given blockSize = {}, and dim() = {}",
344 Dune::BlockVector<Dune::FieldVector<T, blockSize>> returnValue(dim() / blockSize);
345 copyToHost(returnValue);
363 void setZeroAtIndexSet(
const CuVector<int>& indexSet);
366 std::string toDebugString()
368 std::vector<T> v = asStdVector();
369 std::string res =
"";
371 res += std::to_string(element) +
" ";
373 res += std::to_string(v[v.size()-1]);
378 T* m_dataOnDevice =
nullptr;
382 const int m_numberOfElements;
383 detail::CuBlasHandle& m_cuBlasHandle;
385 void assertSameSize(
const CuVector<T>& other)
const;
386 void assertSameSize(
int size)
const;
388 void assertHasElements()
const;