20#ifndef OPM_HYPRE_GPU_TRANSFERS_HPP
21#define OPM_HYPRE_GPU_TRANSFERS_HPP
23#include <opm/simulators/linalg/gpuistl/hypreinterface/HypreDataStructures.hpp>
24#include <opm/simulators/linalg/gpuistl/hypreinterface/HypreErrorHandling.hpp>
28#include <opm/simulators/linalg/gpuistl_hip/GpuSparseMatrixWrapper.hpp>
29#include <opm/simulators/linalg/gpuistl_hip/GpuVector.hpp>
31#include <opm/simulators/linalg/gpuistl/GpuSparseMatrixWrapper.hpp>
32#include <opm/simulators/linalg/gpuistl/GpuVector.hpp>
37#include <_hypre_utilities.h>
42template <
typename VectorType>
43void setContinuousGpuVectorForHypre(
const VectorType&
v,
44 std::vector<HYPRE_Real>& continuous_vector_values,
45 const std::vector<int>& local_hypre_to_local_dune);
47template <
typename VectorType>
48void setGpuVectorFromContinuousVector(VectorType&
v,
49 const std::vector<HYPRE_Real>& continuous_vector_values,
50 const std::vector<int>& local_hypre_to_local_dune);
51#if HYPRE_USING_CUDA || HYPRE_USING_HIP
56template <
typename VectorType>
65 const int N =
static_cast<int>(
host_arrays.indices.size());
66 using T =
typename VectorType::field_type;
77 setContinuousGpuVectorForHypre(
99 setContinuousGpuVectorForHypre(
113template <
typename VectorType>
122 const int N =
static_cast<int>(
host_arrays.indices.size());
123 using T =
typename VectorType::field_type;
141 setGpuVectorFromContinuousVector(
160 setGpuVectorFromContinuousVector(
171template <
typename MatrixType>
181 using T =
typename MatrixType::field_type;
182 const T* values =
gpu_matrix.getNonZeroValues().data();
208template <
typename VectorType>
210setContinuousGpuVectorForHypre(
const VectorType&
v,
211 std::vector<HYPRE_Real>& continuous_vector_values,
212 const std::vector<int>& local_hypre_to_local_dune)
217 for (
size_t i = 0; i < local_hypre_to_local_dune.size(); ++i) {
218 continuous_vector_values[i] =
host_values[local_hypre_to_local_dune[i]];
222template <
typename VectorType>
224setGpuVectorFromContinuousVector(VectorType&
v,
225 const std::vector<HYPRE_Real>& continuous_vector_values,
226 const std::vector<int>& local_hypre_to_local_dune)
230 for (
size_t i = 0; i < local_hypre_to_local_dune.size(); ++i) {
231 host_values[local_hypre_to_local_dune[i]] = continuous_vector_values[i];
Unified interface for Hypre operations with both CPU and GPU data structures.
Definition HypreInterface.hpp:61
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240