GpuVector.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_GPUVECTOR_HEADER_HPP
20#define OPM_GPUVECTOR_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::gpuistl
33{
34
64template <typename T>
65class GpuVector
66{
67public:
68 using field_type = T;
69 using size_type = size_t;
70
71
80 GpuVector(const GpuVector<T>& other);
81
93 explicit GpuVector(const std::vector<T>& data);
94
103 GpuVector& operator=(const GpuVector<T>& other);
104
116 template<int BlockDimension>
117 explicit GpuVector(const Dune::BlockVector<Dune::FieldVector<T, BlockDimension>>& bvector)
118 : GpuVector(bvector.dim())
119 {
120 copyFromHost(bvector);
121 }
122
130 GpuVector& operator=(T scalar);
131
139 explicit GpuVector(const size_t numberOfElements);
140
141
153 GpuVector(const T* dataOnHost, const size_t numberOfElements);
154
158 virtual ~GpuVector();
159
163 T* data();
164
168 const T* data() const;
169
177 template <int BlockDimension>
178 void copyFromHost(const Dune::BlockVector<Dune::FieldVector<T, BlockDimension>>& bvector)
179 {
180 // TODO: [perf] vector.dim() can be replaced by bvector.N() * BlockDimension
181 if (detail::to_size_t(m_numberOfElements) != bvector.dim()) {
182 OPM_THROW(std::runtime_error,
183 fmt::format("Given incompatible vector size. GpuVector has size {}, \n"
184 "however, BlockVector has N() = {}, and dim = {}.",
185 m_numberOfElements,
186 bvector.N(),
187 bvector.dim()));
188 }
189 const auto dataPointer = static_cast<const T*>(&(bvector[0][0]));
190 copyFromHost(dataPointer, m_numberOfElements);
191 }
192
200 template <int BlockDimension>
201 void copyToHost(Dune::BlockVector<Dune::FieldVector<T, BlockDimension>>& bvector) const
202 {
203 // TODO: [perf] vector.dim() can be replaced by bvector.N() * BlockDimension
204 if (detail::to_size_t(m_numberOfElements) != bvector.dim()) {
205 OPM_THROW(std::runtime_error,
206 fmt::format("Given incompatible vector size. GpuVector has size {},\n however, the BlockVector "
207 "has has N() = {}, and dim() = {}.",
208 m_numberOfElements,
209 bvector.N(),
210 bvector.dim()));
211 }
212 const auto dataPointer = static_cast<T*>(&(bvector[0][0]));
213 copyToHost(dataPointer, m_numberOfElements);
214 }
215
223 void copyFromHost(const T* dataPointer, size_t numberOfElements);
224 void copyFromHost(const T* dataPointer, size_t numberOfElements, cudaStream_t stream);
225
233 void copyToHost(T* dataPointer, size_t numberOfElements) const;
234
242 void copyFromHost(const std::vector<T>& data);
243
251 void copyToHost(std::vector<T>& data) const;
252
257 void copyFromDeviceToDevice(const GpuVector<T>& other) const;
258
259 void prepareSendBuf(GpuVector<T>& buffer, const GpuVector<int>& indexSet) const;
260 void syncFromRecvBuf(GpuVector<T>& buffer, const GpuVector<int>& indexSet) const;
261
270 GpuVector<T>& operator*=(const T& scalar);
271
280 GpuVector<T>& axpy(T alpha, const GpuVector<T>& y);
281
288 GpuVector<T>& operator+=(const GpuVector<T>& other);
289
296 GpuVector<T>& operator-=(const GpuVector<T>& other);
297
306 T dot(const GpuVector<T>& other) const;
307
315 T two_norm() const;
316
322 T dot(const GpuVector<T>& other, const GpuVector<int>& indexSet, GpuVector<T>& buffer) const;
323
329 T two_norm(const GpuVector<int>& indexSet, GpuVector<T>& buffer) const;
330
331
337 T dot(const GpuVector<T>& other, const GpuVector<int>& indexSet) const;
338
344 T two_norm(const GpuVector<int>& indexSet) const;
345
346
351 size_type dim() const;
352
353
358 std::vector<T> asStdVector() const;
359
364 template <int blockSize>
365 Dune::BlockVector<Dune::FieldVector<T, blockSize>> asDuneBlockVector() const
366 {
367 OPM_ERROR_IF(dim() % blockSize != 0,
368 fmt::format("blockSize is not a multiple of dim(). Given blockSize = {}, and dim() = {}",
369 blockSize,
370 dim()));
371
372 Dune::BlockVector<Dune::FieldVector<T, blockSize>> returnValue(dim() / blockSize);
373 copyToHost(returnValue);
374 return returnValue;
375 }
376
377
391 void setZeroAtIndexSet(const GpuVector<int>& indexSet);
392
393 // Slow method that creates a string representation of a GpuVector for debug purposes
394 std::string toDebugString()
395 {
396 std::vector<T> v = asStdVector();
397 std::string res = "";
398 for (T element : v){
399 res += std::to_string(element) + " ";
400 }
401 res += std::to_string(v[v.size()-1]);
402 return res;
403 }
404
405private:
406 T* m_dataOnDevice = nullptr;
407
408 // Note that we store this as int to make sure we are always cublas compatible.
409 // This gives the added benefit that a size_t to int conversion error occurs during construction.
410 const int m_numberOfElements;
411 detail::CuBlasHandle& m_cuBlasHandle;
412
413 void assertSameSize(const GpuVector<T>& other) const;
414 void assertSameSize(int size) const;
415
416 void assertHasElements() const;
417};
418
419} // namespace Opm::gpuistl
420#endif
static constexpr int dim
Definition: structuredgridvanguard.hh:68
void syncFromRecvBuf(T *deviceA, T *buffer, size_t numberOfElements, const int *indices)
void prepareSendBuf(const T *deviceA, T *buffer, size_t numberOfElements, const int *indices)
__host__ __device__ 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:86
Definition: autotuner.hpp:30
std::string toDebugString()
Definition: GpuVector.hpp:394
void setZeroAtIndexSet(const GpuVector< int > &indexSet)
The GpuVector class is a simple (arithmetic) vector class for the GPU.
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)