CuVector.hpp
Go to the documentation of this file.
1/*
2 Copyright 2022-2023 SINTEF AS
3
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18*/
19#ifndef OPM_CUVECTOR_HEADER_HPP
20#define OPM_CUVECTOR_HEADER_HPP
21#include <dune/common/fvector.hh>
22#include <dune/istl/bvector.hh>
23#include <exception>
24#include <fmt/core.h>
25#include <opm/common/ErrorMacros.hpp>
28#include <vector>
29#include <string>
30
31
32namespace Opm::cuistl
33{
34
64template <typename T>
65class CuVector
66{
67public:
68 using field_type = T;
69 using size_type = size_t;
70
71
80 CuVector(const CuVector<T>& other);
81
93 explicit CuVector(const std::vector<T>& data);
94
103 CuVector& operator=(const CuVector<T>& other);
104
112 CuVector& operator=(T scalar);
113
121 explicit CuVector(const size_t numberOfElements);
122
123
135 CuVector(const T* dataOnHost, const size_t numberOfElements);
136
140 virtual ~CuVector();
141
145 T* data();
146
150 const T* data() const;
151
159 template <int BlockDimension>
160 void copyFromHost(const Dune::BlockVector<Dune::FieldVector<T, BlockDimension>>& bvector)
161 {
162 // TODO: [perf] vector.dim() can be replaced by bvector.N() * BlockDimension
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 = {}.",
167 m_numberOfElements,
168 bvector.N(),
169 bvector.dim()));
170 }
171 const auto dataPointer = static_cast<const T*>(&(bvector[0][0]));
172 copyFromHost(dataPointer, m_numberOfElements);
173 }
174
182 template <int BlockDimension>
183 void copyToHost(Dune::BlockVector<Dune::FieldVector<T, BlockDimension>>& bvector) const
184 {
185 // TODO: [perf] vector.dim() can be replaced by bvector.N() * BlockDimension
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() = {}.",
190 m_numberOfElements,
191 bvector.N(),
192 bvector.dim()));
193 }
194 const auto dataPointer = static_cast<T*>(&(bvector[0][0]));
195 copyToHost(dataPointer, m_numberOfElements);
196 }
197
205 void copyFromHost(const T* dataPointer, size_t numberOfElements);
206
214 void copyToHost(T* dataPointer, size_t numberOfElements) const;
215
223 void copyFromHost(const std::vector<T>& data);
224
232 void copyToHost(std::vector<T>& data) const;
233
234 void prepareSendBuf(CuVector<T>& buffer, const CuVector<int>& indexSet) const;
235 void syncFromRecvBuf(CuVector<T>& buffer, const CuVector<int>& indexSet) const;
236
245 CuVector<T>& operator*=(const T& scalar);
246
255 CuVector<T>& axpy(T alpha, const CuVector<T>& y);
256
263 CuVector<T>& operator+=(const CuVector<T>& other);
264
271 CuVector<T>& operator-=(const CuVector<T>& other);
272
281 T dot(const CuVector<T>& other) const;
282
290 T two_norm() const;
291
297 T dot(const CuVector<T>& other, const CuVector<int>& indexSet, CuVector<T>& buffer) const;
298
304 T two_norm(const CuVector<int>& indexSet, CuVector<T>& buffer) const;
305
306
312 T dot(const CuVector<T>& other, const CuVector<int>& indexSet) const;
313
319 T two_norm(const CuVector<int>& indexSet) const;
320
321
326 size_type dim() const;
327
328
333 std::vector<T> asStdVector() const;
334
339 template <int blockSize>
340 Dune::BlockVector<Dune::FieldVector<T, blockSize>> asDuneBlockVector() const
341 {
342 OPM_ERROR_IF(dim() % blockSize != 0,
343 fmt::format("blockSize is not a multiple of dim(). Given blockSize = {}, and dim() = {}",
344 blockSize,
345 dim()));
346
347 Dune::BlockVector<Dune::FieldVector<T, blockSize>> returnValue(dim() / blockSize);
348 copyToHost(returnValue);
349 return returnValue;
350 }
351
352
366 void setZeroAtIndexSet(const CuVector<int>& indexSet);
367
368 // Slow method that creates a string representation of a CuVector for debug purposes
369 std::string toDebugString()
370 {
371 std::vector<T> v = asStdVector();
372 std::string res = "";
373 for (T element : v){
374 res += std::to_string(element) + " ";
375 }
376 res += std::to_string(v[v.size()-1]);
377 return res;
378 }
379
380private:
381 T* m_dataOnDevice = nullptr;
382
383 // Note that we store this as int to make sure we are always cublas compatible.
384 // This gives the added benefit that a size_t to int conversion error occurs during construction.
385 const int m_numberOfElements;
386 detail::CuBlasHandle& m_cuBlasHandle;
387
388 void assertSameSize(const CuVector<T>& other) const;
389 void assertSameSize(int size) const;
390
391 void assertHasElements() const;
392};
393
394} // namespace Opm::cuistl
395#endif
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)