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>
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)
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
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;
234 void prepareSendBuf(CuVector<T>& buffer,
const CuVector<int>& indexSet)
const;
235 void syncFromRecvBuf(CuVector<T>& buffer,
const CuVector<int>& indexSet)
const;
245 CuVector<T>& operator*=(
const T& scalar);
255 CuVector<T>& axpy(T alpha,
const CuVector<T>& y);
263 CuVector<T>& operator+=(
const CuVector<T>& other);
271 CuVector<T>& operator-=(
const CuVector<T>& other);
281 T dot(
const CuVector<T>& other)
const;
297 T dot(
const CuVector<T>& other,
const CuVector<int>& indexSet, CuVector<T>& buffer)
const;
304 T two_norm(
const CuVector<int>& indexSet, CuVector<T>& buffer)
const;
312 T dot(
const CuVector<T>& other,
const CuVector<int>& indexSet)
const;
319 T two_norm(
const CuVector<int>& indexSet)
const;
326 size_type dim()
const;
333 std::vector<T> asStdVector()
const;
339 template <
int blockSize>
340 Dune::BlockVector<Dune::FieldVector<T, blockSize>> asDuneBlockVector()
const
342 OPM_ERROR_IF(dim() % blockSize != 0,
343 fmt::format(
"blockSize is not a multiple of dim(). Given blockSize = {}, and dim() = {}",
347 Dune::BlockVector<Dune::FieldVector<T, blockSize>> returnValue(dim() / blockSize);
348 copyToHost(returnValue);
371 std::vector<T> v = asStdVector();
372 std::string res =
"";
381 T* m_dataOnDevice =
nullptr;
385 const int m_numberOfElements;
386 detail::CuBlasHandle& m_cuBlasHandle;
388 void assertSameSize(
const CuVector<T>& other)
const;
389 void assertSameSize(
int size)
const;
391 void assertHasElements()
const;
void prepareSendBuf(const T *deviceA, T *buffer, size_t numberOfElements, const int *indices)
void syncFromRecvBuf(T *deviceA, T *buffer, size_t numberOfElements, const int *indices)
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:85
Definition: CuBlockPreconditioner.hpp:29
void setZeroAtIndexSet(const CuVector< int > &indexSet)
The CuVector class is a simple (arithmetic) vector class for the GPU.
std::string toDebugString()
Definition: CuVector.hpp:369
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)