20#ifndef OPM_HYPRE_SETUP_HPP
21#define OPM_HYPRE_SETUP_HPP
23#include <opm/simulators/linalg/PropertyTree.hpp>
24#include <opm/simulators/linalg/gpuistl/detail/gpu_type_detection.hpp>
25#include <opm/simulators/linalg/gpuistl/hypreinterface/HypreDataStructures.hpp>
26#include <opm/simulators/linalg/gpuistl/hypreinterface/HypreErrorHandling.hpp>
28#include <dune/istl/owneroverlapcopy.hh>
29#include <dune/istl/paamg/graph.hh>
30#include <dune/istl/paamg/pinfo.hh>
31#include <dune/istl/repartition.hh>
35#include <opm/simulators/linalg/gpuistl_hip/GpuSparseMatrixWrapper.hpp>
37#include <opm/simulators/linalg/gpuistl/GpuSparseMatrixWrapper.hpp>
42#include <HYPRE_parcsr_ls.h>
43#include <_hypre_utilities.h>
53template <
typename T,
bool ForceLegacy>
54SparsityPattern setupSparsityPatternFromGpuMatrix(
const GpuSparseMatrixWrapper<T, ForceLegacy>& gpu_matrix,
55 const ParallelInfo& par_info,
58template <
typename T,
bool ForceLegacy>
59std::vector<HYPRE_Int> computeRowIndexesWithMappingGpu(
const GpuSparseMatrixWrapper<T, ForceLegacy>& gpu_matrix,
60 const std::vector<int>& local_dune_to_local_hypre);
66template <
typename CommType,
typename MatrixType>
69template <
typename MatrixType>
73template <
typename MatrixType>
74std::vector<HYPRE_Int> computeRowIndexesWithMappingCpu(
const MatrixType& matrix,
75 const std::vector<HYPRE_Int>& ncols,
76 const std::vector<int>& local_dune_to_local_hypre,
79template <
typename MatrixType>
80std::vector<HYPRE_Int> computeRowIndexesWithMappingCpu(
const MatrixType& matrix,
81 const std::vector<int>& local_dune_to_local_hypre);
92#if HYPRE_USING_CUDA || HYPRE_USING_HIP
93 if (use_gpu_backend) {
94 OPM_HYPRE_SAFE_CALL(HYPRE_SetMemoryLocation(HYPRE_MEMORY_DEVICE));
95 OPM_HYPRE_SAFE_CALL(HYPRE_SetExecutionPolicy(HYPRE_EXEC_DEVICE));
97 OPM_HYPRE_SAFE_CALL(HYPRE_SetSpGemmUseVendor(
false));
99 OPM_HYPRE_SAFE_CALL(HYPRE_SetUseGpuRand(1));
100 OPM_HYPRE_SAFE_CALL(HYPRE_DeviceInitialize());
101 OPM_HYPRE_SAFE_CALL(HYPRE_PrintDeviceInfo());
105 OPM_HYPRE_SAFE_CALL(HYPRE_SetMemoryLocation(HYPRE_MEMORY_HOST));
106 OPM_HYPRE_SAFE_CALL(HYPRE_SetExecutionPolicy(HYPRE_EXEC_HOST));
120 OPM_HYPRE_SAFE_CALL(HYPRE_BoomerAMGCreate(&solver));
136 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetPrintLevel(solver, prm.
get<
int>(
"print_level", 0)));
137 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetMaxIter(solver, prm.
get<
int>(
"max_iter", 1)));
138 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetStrongThreshold(solver, prm.
get<
double>(
"strong_threshold", 0.5)));
139 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetAggTruncFactor(solver, prm.
get<
double>(
"agg_trunc_factor", 0.3)));
140 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetInterpType(solver, prm.
get<
int>(
"interp_type", 6)));
141 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetMaxLevels(solver, prm.
get<
int>(
"max_levels", 15)));
142 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetTol(solver, prm.
get<
double>(
"tolerance", 0.0)));
144 if (use_gpu_backend) {
145 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetRelaxType(solver, prm.
get<
int>(
"relax_type", 16)));
146 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetCoarsenType(solver, prm.
get<
int>(
"coarsen_type", 8)));
147 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetAggNumLevels(solver, prm.
get<
int>(
"agg_num_levels", 0)));
148 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetAggInterpType(solver, prm.
get<
int>(
"agg_interp_type", 6)));
150 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetKeepTranspose(solver,
true));
152 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetRelaxType(solver, prm.
get<
int>(
"relax_type", 13)));
153 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetCoarsenType(solver, prm.
get<
int>(
"coarsen_type", 10)));
154 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetAggNumLevels(solver, prm.
get<
int>(
"agg_num_levels", 1)));
155 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetAggInterpType(solver, prm.
get<
int>(
"agg_interp_type", 4)));
168template <
typename CommType>
172 HYPRE_IJMatrix matrix;
174 if constexpr (std::is_same_v<CommType, Dune::Amg::SequentialInformation>) {
175 mpi_comm = MPI_COMM_SELF;
177 mpi_comm = comm.communicator();
180 HYPRE_IJMatrixCreate(mpi_comm, dof_offset, dof_offset + (N - 1), dof_offset, dof_offset + (N - 1), &matrix));
181 OPM_HYPRE_SAFE_CALL(HYPRE_IJMatrixSetObjectType(matrix, HYPRE_PARCSR));
182 OPM_HYPRE_SAFE_CALL(HYPRE_IJMatrixInitialize(matrix));
195template <
typename CommType>
199 HYPRE_IJVector vector;
201 if constexpr (std::is_same_v<CommType, Dune::Amg::SequentialInformation>) {
202 mpi_comm = MPI_COMM_SELF;
204 mpi_comm = comm.communicator();
206 OPM_HYPRE_SAFE_CALL(HYPRE_IJVectorCreate(mpi_comm, dof_offset, dof_offset + (N - 1), &vector));
207 OPM_HYPRE_SAFE_CALL(HYPRE_IJVectorSetObjectType(vector, HYPRE_PARCSR));
208 OPM_HYPRE_SAFE_CALL(HYPRE_IJVectorInitialize(vector));
222 OPM_HYPRE_SAFE_CALL(HYPRE_BoomerAMGDestroy(solver));
236 OPM_HYPRE_SAFE_CALL(HYPRE_IJMatrixDestroy(matrix));
250 OPM_HYPRE_SAFE_CALL(HYPRE_IJVectorDestroy(vector));
261template <
typename CommType,
typename MatrixType>
265 if constexpr (std::is_same_v<CommType, Dune::Amg::SequentialInformation>) {
361template <
typename CommType,
typename MatrixType>
366 const auto& collective_comm = comm.communicator();
373 if (!(matrix.N() == comm.indexSet().size())) {
376 const_cast<CommType&
>(comm).buildGlobalLookup(matrix.N());
377 Dune::Amg::MatrixGraph<MatrixType> graph(
const_cast<MatrixType&
>(matrix));
378 Dune::fillIndexSetHoles(graph,
const_cast<CommType&
>(comm));
379 assert(matrix.N() == comm.indexSet().size());
385 for (
const auto& ind : comm.indexSet()) {
386 int local_ind = ind.local().local();
387 if (ind.local().attribute() == Dune::OwnerOverlapCopyAttributeSet::owner) {
399 bool owner_first =
true;
400 bool visited_ghost =
false;
401 std::size_t count = 0;
405 visited_ghost =
true;
413 owner_first = owner_first && !visited_ghost;
425 std::vector<int> dof_counts_per_process(collective_comm.size());
426 collective_comm.allgather(&info.
N_owned, 1, dof_counts_per_process.data());
429 info.
dof_offset = std::accumulate(dof_counts_per_process.begin(),
430 dof_counts_per_process.begin() + collective_comm.rank(),
444 if (collective_comm.rank() > 0) {
463template <
typename MatrixType>
467#if HYPRE_USING_CUDA || HYPRE_USING_HIP
469 return setupSparsityPatternFromGpuMatrix(matrix, par_info, owner_first);
485template <
typename MatrixType>
489 SparsityPattern pattern;
493 std::size_t cols_size = 0;
495 for (
auto row = matrix.begin(); row != matrix.end(); ++row) {
496 const int rowIdx = row.index();
498 cols_size += row->size();
501 pattern.
nnz = cols_size;
504 pattern.
nnz = matrix.nonzeroes();
510 pattern.
cols.resize(pattern.
nnz);
513 for (
auto row = matrix.begin(); row != matrix.end(); ++row) {
514 const int rind = row.index();
519 if (owner_first && local_rowIdx < 0) {
523 if (local_rowIdx >= 0) {
526 pattern.
rows[local_rowIdx] = global_rowIdx;
527 pattern.
ncols[local_rowIdx] = row->size();
531 for (
auto col = row->begin(); col != row->end(); ++col) {
533 assert(global_colIdx >= 0);
534 pattern.
cols[pos++] = global_colIdx;
541#if HYPRE_USING_CUDA || HYPRE_USING_HIP
550template <
typename T,
bool ForceLegacy>
560 std::size_t cols_size = 0;
562 auto host_row_ptrs = gpu_matrix.
getRowIndices().asStdVector();
563 for (
int rind = 0; rind < static_cast<int>(gpu_matrix.
N()); ++rind) {
565 const int row_start = host_row_ptrs[rind];
566 const int row_end = host_row_ptrs[rind + 1];
567 cols_size += (row_end - row_start);
570 pattern.
nnz = cols_size;
579 pattern.
cols.resize(pattern.
nnz);
582 auto host_row_ptrs = gpu_matrix.
getRowIndices().asStdVector();
586 for (
int rind = 0; rind < static_cast<int>(gpu_matrix.
N()); ++rind) {
591 if (owner_first && local_rowIdx < 0) {
595 const int row_start = host_row_ptrs[rind];
596 const int row_end = host_row_ptrs[rind + 1];
597 const int num_cols = row_end - row_start;
599 if (local_rowIdx >= 0) {
602 pattern.
rows[local_rowIdx] = global_rowIdx;
603 pattern.
ncols[local_rowIdx] = num_cols;
607 for (
int col_idx = row_start; col_idx < row_end; ++col_idx) {
608 const int colIdx = host_col_indices[col_idx];
610 assert(global_colIdx >= 0);
611 pattern.
cols[pos++] = global_colIdx;
632template <
typename MatrixType>
633std::vector<HYPRE_Int>
635 const std::vector<HYPRE_Int>& ncols,
636 const std::vector<int>& local_dune_to_local_hypre,
641 std::vector<HYPRE_Int> row_indexes(ncols.size());
643 for (std::size_t i = 1; i < ncols.size(); ++i) {
644 row_indexes[i] = row_indexes[i - 1] + ncols[i - 1];
649#if HYPRE_USING_CUDA || HYPRE_USING_HIP
651 return computeRowIndexesWithMappingGpu(matrix, local_dune_to_local_hypre);
655 return computeRowIndexesWithMappingCpu(matrix, local_dune_to_local_hypre);
679template <
typename MatrixType>
680std::vector<HYPRE_Int>
681computeRowIndexesWithMappingCpu(
const MatrixType& matrix,
const std::vector<int>& local_dune_to_local_hypre)
683 const int N = std::count_if(
684 local_dune_to_local_hypre.begin(), local_dune_to_local_hypre.end(), [](
int val) { return val >= 0; });
685 std::vector<HYPRE_Int> row_indexes(N);
686 int data_position = 0;
688 for (
auto row = matrix.begin(); row != matrix.end(); ++row) {
689 const int dune_row_idx = row.index();
690 const int hypre_row_idx = local_dune_to_local_hypre[dune_row_idx];
692 if (hypre_row_idx >= 0) {
695 row_indexes[hypre_row_idx] = data_position;
698 data_position += row->size();
703#if HYPRE_USING_CUDA || HYPRE_USING_HIP
715template <
typename T,
bool ForceLegacy>
716std::vector<HYPRE_Int>
718 const std::vector<int>& local_dune_to_local_hypre)
720 const int N = std::count_if(
721 local_dune_to_local_hypre.begin(), local_dune_to_local_hypre.end(), [](
int val) { return val >= 0; });
722 std::vector<HYPRE_Int> row_indexes(N);
725 auto host_row_ptrs = gpu_matrix.
getRowIndices().asStdVector();
728 for (
int dune_row_idx = 0; dune_row_idx < static_cast<int>(gpu_matrix.
N()); ++dune_row_idx) {
729 const int hypre_row_idx = local_dune_to_local_hypre[dune_row_idx];
731 if (hypre_row_idx >= 0) {
734 row_indexes[hypre_row_idx] = host_row_ptrs[dune_row_idx];
Hierarchical collection of key/value pairs.
Definition PropertyTree.hpp:39
T get(const std::string &key) const
Retrieve property value given hierarchical property key.
Definition PropertyTree.cpp:59
The GpuSparseMatrixWrapper Checks CUDA/HIP version and dispatches a version either using the old or t...
Definition GpuSparseMatrixWrapper.hpp:62
GpuVector< int > & getRowIndices()
getRowIndices returns the row indices used to represent the BSR structure.
Definition GpuSparseMatrixWrapper.hpp:271
std::size_t N() const
N returns the number of rows (which is equal to the number of columns).
Definition GpuSparseMatrixWrapper.hpp:232
std::size_t nonzeroes() const
nonzeroes behaves as the Dune::BCRSMatrix::nonzeros() function and returns the number of non zero blo...
Definition GpuSparseMatrixWrapper.hpp:241
GpuVector< int > & getColumnIndices()
getColumnIndices returns the column indices used to represent the BSR structure.
Definition GpuSparseMatrixWrapper.hpp:291
Unified interface for Hypre operations with both CPU and GPU data structures.
Definition HypreInterface.hpp:61
void destroySolver(HYPRE_Solver solver)
Destroy Hypre solver.
Definition HypreSetup.hpp:219
SparsityPattern setupSparsityPatternFromCpuMatrix(const MatrixType &matrix, const ParallelInfo &par_info, bool owner_first)
Setup sparsity pattern from CPU matrix (BCRSMatrix).
Definition HypreSetup.hpp:487
ParallelInfo setupHypreParallelInfo(const CommType &comm, const MatrixType &matrix)
Setup parallel information for Hypre (automatically detects serial/parallel).
Definition HypreSetup.hpp:263
void destroyMatrix(HYPRE_IJMatrix matrix)
Destroy Hypre matrix.
Definition HypreSetup.hpp:233
HYPRE_Solver createAMGSolver()
Create Hypre solver (BoomerAMG).
Definition HypreSetup.hpp:117
std::vector< HYPRE_Int > computeRowIndexes(const MatrixType &matrix, const std::vector< HYPRE_Int > &ncols, const std::vector< int > &local_dune_to_local_hypre, bool owner_first)
Compute row indexes for HYPRE_IJMatrixSetValues2.
Definition HypreSetup.hpp:634
HYPRE_IJMatrix createMatrix(HYPRE_Int N, HYPRE_Int dof_offset, const CommType &comm)
Create Hypre matrix.
Definition HypreSetup.hpp:170
ParallelInfo setupHypreParallelInfoParallel(const CommType &comm, const MatrixType &matrix)
Create mappings between Dune and HYPRE indexing for parallel decomposition.
Definition HypreSetup.hpp:363
SparsityPattern setupSparsityPattern(const MatrixType &matrix, const ParallelInfo &par_info, bool owner_first)
Setup sparsity pattern from matrix (automatically detects CPU/GPU type).
Definition HypreSetup.hpp:465
void initialize(bool use_gpu_backend)
Initialize the Hypre library and set memory/execution policy.
Definition HypreSetup.hpp:89
HYPRE_IJVector createVector(HYPRE_Int N, HYPRE_Int dof_offset, const CommType &comm)
Create Hypre vector.
Definition HypreSetup.hpp:197
void destroyVector(HYPRE_IJVector vector)
Destroy Hypre vector.
Definition HypreSetup.hpp:247
void setSolverParameters(HYPRE_Solver solver, const PropertyTree &prm, bool use_gpu_backend)
Set solver parameters from property tree.
Definition HypreSetup.hpp:133
ParallelInfo setupHypreParallelInfoSerial(HYPRE_Int N)
Setup parallel information for Hypre in serial case.
Definition HypreSetup.hpp:279
Parallel domain decomposition information for HYPRE-Dune interface.
Definition HypreDataStructures.hpp:37
bool owner_first
Whether owned DOFs appear first in local Dune ordering.
Definition HypreDataStructures.hpp:77
std::vector< int > local_dune_to_global_hypre
Mapping from local Dune indices to global HYPRE indices.
Definition HypreDataStructures.hpp:51
std::vector< int > local_dune_to_local_hypre
Mapping from local Dune indices to local HYPRE indices.
Definition HypreDataStructures.hpp:44
std::vector< int > local_hypre_to_local_dune
Mapping from local HYPRE indices to local Dune indices.
Definition HypreDataStructures.hpp:59
HYPRE_Int dof_offset
Global index offset for this process's owned DOFs.
Definition HypreDataStructures.hpp:69
HYPRE_Int N_owned
Number of DOFs owned by this MPI process.
Definition HypreDataStructures.hpp:62
Compressed Sparse Row (CSR) sparsity pattern for HYPRE matrix assembly.
Definition HypreDataStructures.hpp:86
std::vector< HYPRE_BigInt > rows
Global row indices for owned rows (size: N_owned).
Definition HypreDataStructures.hpp:91
std::vector< HYPRE_Int > ncols
Non-zero entries per owned row (size: N_owned).
Definition HypreDataStructures.hpp:88
HYPRE_Int nnz
Number of non-zero entries in matrix.
Definition HypreDataStructures.hpp:97
std::vector< HYPRE_BigInt > cols
Global column indices in CSR format (size: nnz).
Definition HypreDataStructures.hpp:94
Type trait to detect if a type is a GPU type.
Definition gpu_type_detection.hpp:40