hypreinterface/HypreSetup.hpp
Go to the documentation of this file.
1/*
2 Copyright 2025 Equinor ASA
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
20#ifndef OPM_HYPRE_SETUP_HPP
21#define OPM_HYPRE_SETUP_HPP
22
26
27#include <dune/istl/owneroverlapcopy.hh>
28#include <dune/istl/paamg/graph.hh>
29#include <dune/istl/paamg/pinfo.hh>
30#include <dune/istl/repartition.hh>
31
32#if HYPRE_USING_CUDA
36#elif HYPRE_USING_HIP
37#include <opm/simulators/linalg/gpuistl_hip/detail/gpu_type_detection.hpp>
38#include <opm/simulators/linalg/gpuistl_hip/GpuSparseMatrixWrapper.hpp>
39#include <opm/simulators/linalg/gpuistl_hip/hypreinterface/HypreSetup.hpp>
40#endif
41
42#include <HYPRE.h>
43#include <HYPRE_parcsr_ls.h>
44#include <_hypre_utilities.h>
45
46#include <algorithm>
47#include <cstddef>
48#include <numeric>
49
51{
52
53// Serial helper functions
54ParallelInfo setupHypreParallelInfoSerial(HYPRE_Int N);
55
56// Parallel helper functions
57template <typename CommType, typename MatrixType>
58ParallelInfo
59setupHypreParallelInfoParallel(const CommType& comm, const MatrixType& matrix);
60
61template <typename MatrixType>
62SparsityPattern
63setupSparsityPatternFromCpuMatrix(const MatrixType& matrix,
64 const ParallelInfo& par_info,
65 bool owner_first);
66
67template <typename MatrixType>
68std::vector<HYPRE_Int> computeRowIndexesWithMappingCpu(const MatrixType& matrix,
69 const std::vector<HYPRE_Int>& ncols,
70 const std::vector<int>& local_dune_to_local_hypre,
71 bool owner_first);
72
73template <typename MatrixType>
74std::vector<HYPRE_Int> computeRowIndexesWithMappingCpu(const MatrixType& matrix,
75 const std::vector<int>& local_dune_to_local_hypre);
76
83inline void
84initialize([[maybe_unused]] bool use_gpu_backend)
85{
86 // Set memory location and execution policy
87#if HYPRE_USING_CUDA || HYPRE_USING_HIP
88 if (use_gpu_backend) {
89 OPM_HYPRE_SAFE_CALL(HYPRE_SetMemoryLocation(HYPRE_MEMORY_DEVICE));
90 OPM_HYPRE_SAFE_CALL(HYPRE_SetExecutionPolicy(HYPRE_EXEC_DEVICE));
91 // use hypre's SpGEMM instead of vendor implementation
92 OPM_HYPRE_SAFE_CALL(HYPRE_SetSpGemmUseVendor(false));
93 // use cuRand for PMIS
94 OPM_HYPRE_SAFE_CALL(HYPRE_SetUseGpuRand(1));
95 OPM_HYPRE_SAFE_CALL(HYPRE_DeviceInitialize());
96 OPM_HYPRE_SAFE_CALL(HYPRE_PrintDeviceInfo());
97 } else
98#endif
99 {
100 OPM_HYPRE_SAFE_CALL(HYPRE_SetMemoryLocation(HYPRE_MEMORY_HOST));
101 OPM_HYPRE_SAFE_CALL(HYPRE_SetExecutionPolicy(HYPRE_EXEC_HOST));
102 }
103}
104
111inline HYPRE_Solver
113{
114 HYPRE_Solver solver;
115 OPM_HYPRE_SAFE_CALL(HYPRE_BoomerAMGCreate(&solver));
116 return solver;
117}
118
127inline void
128setSolverParameters(HYPRE_Solver solver, const PropertyTree& prm, bool use_gpu_backend)
129{
130 // Set parameters from property tree with defaults
131 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetPrintLevel(solver, prm.get<int>("print_level", 0)));
132 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetMaxIter(solver, prm.get<int>("max_iter", 1)));
133 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetStrongThreshold(solver, prm.get<double>("strong_threshold", 0.5)));
134 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetAggTruncFactor(solver, prm.get<double>("agg_trunc_factor", 0.3)));
135 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetInterpType(solver, prm.get<int>("interp_type", 6)));
136 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetMaxLevels(solver, prm.get<int>("max_levels", 15)));
137 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetTol(solver, prm.get<double>("tolerance", 0.0)));
138
139 if (use_gpu_backend) {
140 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetRelaxType(solver, prm.get<int>("relax_type", 16)));
141 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetCoarsenType(solver, prm.get<int>("coarsen_type", 8)));
142 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetAggNumLevels(solver, prm.get<int>("agg_num_levels", 0)));
143 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetAggInterpType(solver, prm.get<int>("agg_interp_type", 6)));
144 // Keep transpose to avoid SpMTV
145 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetKeepTranspose(solver, true));
146 } else {
147 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetRelaxType(solver, prm.get<int>("relax_type", 13)));
148 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetCoarsenType(solver, prm.get<int>("coarsen_type", 10)));
149 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetAggNumLevels(solver, prm.get<int>("agg_num_levels", 1)));
150 HYPRE_SAFE_CALL(HYPRE_BoomerAMGSetAggInterpType(solver, prm.get<int>("agg_interp_type", 4)));
151 }
152}
153
163template <typename CommType>
164HYPRE_IJMatrix
165createMatrix(HYPRE_Int N, HYPRE_Int dof_offset, const CommType& comm)
166{
167 HYPRE_IJMatrix matrix;
168 MPI_Comm mpi_comm;
169 if constexpr (std::is_same_v<CommType, Dune::Amg::SequentialInformation>) {
170 mpi_comm = MPI_COMM_SELF;
171 } else {
172 mpi_comm = comm.communicator();
173 }
175 HYPRE_IJMatrixCreate(mpi_comm, dof_offset, dof_offset + (N - 1), dof_offset, dof_offset + (N - 1), &matrix));
176 OPM_HYPRE_SAFE_CALL(HYPRE_IJMatrixSetObjectType(matrix, HYPRE_PARCSR));
177 OPM_HYPRE_SAFE_CALL(HYPRE_IJMatrixInitialize(matrix));
178 return matrix;
179}
180
190template <typename CommType>
191HYPRE_IJVector
192createVector(HYPRE_Int N, HYPRE_Int dof_offset, const CommType& comm)
193{
194 HYPRE_IJVector vector;
195 MPI_Comm mpi_comm;
196 if constexpr (std::is_same_v<CommType, Dune::Amg::SequentialInformation>) {
197 mpi_comm = MPI_COMM_SELF;
198 } else {
199 mpi_comm = comm.communicator();
200 }
201 OPM_HYPRE_SAFE_CALL(HYPRE_IJVectorCreate(mpi_comm, dof_offset, dof_offset + (N - 1), &vector));
202 OPM_HYPRE_SAFE_CALL(HYPRE_IJVectorSetObjectType(vector, HYPRE_PARCSR));
203 OPM_HYPRE_SAFE_CALL(HYPRE_IJVectorInitialize(vector));
204 return vector;
205}
206
213inline void
214destroySolver(HYPRE_Solver solver)
215{
216 if (solver) {
217 OPM_HYPRE_SAFE_CALL(HYPRE_BoomerAMGDestroy(solver));
218 }
219}
220
227inline void
228destroyMatrix(HYPRE_IJMatrix matrix)
229{
230 if (matrix) {
231 OPM_HYPRE_SAFE_CALL(HYPRE_IJMatrixDestroy(matrix));
232 }
233}
234
241inline void
242destroyVector(HYPRE_IJVector vector)
243{
244 if (vector) {
245 OPM_HYPRE_SAFE_CALL(HYPRE_IJVectorDestroy(vector));
246 }
247}
248
256template <typename CommType, typename MatrixType>
257ParallelInfo
258setupHypreParallelInfo(const CommType& comm, const MatrixType& matrix)
259{
260 if constexpr (std::is_same_v<CommType, Dune::Amg::SequentialInformation>) {
261 return setupHypreParallelInfoSerial(static_cast<HYPRE_Int>(matrix.N()));
262 } else {
263 return setupHypreParallelInfoParallel(comm, matrix);
264 }
265}
266
273inline ParallelInfo
275{
276 ParallelInfo info;
277 info.N_owned = N;
278
279 info.local_dune_to_local_hypre.resize(N);
280 info.local_dune_to_global_hypre.resize(N);
281 info.local_hypre_to_local_dune.resize(N);
282
283 std::iota(info.local_dune_to_local_hypre.begin(), info.local_dune_to_local_hypre.end(), 0);
284 std::iota(info.local_hypre_to_local_dune.begin(), info.local_hypre_to_local_dune.end(), 0);
285 std::iota(info.local_dune_to_global_hypre.begin(), info.local_dune_to_global_hypre.end(), 0);
286
287 info.dof_offset = 0;
288 info.owner_first = true;
289
290 return info;
291}
292
356template <typename CommType, typename MatrixType>
357ParallelInfo
358setupHypreParallelInfoParallel(const CommType& comm, const MatrixType& matrix)
359{
360 ParallelInfo info;
361 const auto& collective_comm = comm.communicator();
362
363 // Initialize mapping arrays to not owned (-1) state
364 info.local_dune_to_local_hypre.resize(comm.indexSet().size(), -1);
365 info.local_dune_to_global_hypre.resize(comm.indexSet().size(), -1);
366
367 // Handle edge case: ensure index set covers all matrix rows
368 if (!(matrix.N() == comm.indexSet().size())) {
369 // in OPM this will likely not be trigged
370 // ensure no holes in index sett
371 const_cast<CommType&>(comm).buildGlobalLookup(matrix.N()); // need?
372 Dune::Amg::MatrixGraph<MatrixType> graph(const_cast<MatrixType&>(matrix)); // do not know why not const ref is sufficient
373 Dune::fillIndexSetHoles(graph, const_cast<CommType&>(comm));
374 assert(matrix.N() == comm.indexSet().size());
375 }
376
377 // STEP 1: Ownership Detection
378 // Scan Dune's index set to identify which DOFs this process owns
379 // Note: iteration order in index set is NOT sequential by local index
380 for (const auto& ind : comm.indexSet()) {
381 int local_ind = ind.local().local();
382 if (ind.local().attribute() == Dune::OwnerOverlapCopyAttributeSet::owner) {
383 // Mark as owned (temporarily use 1, will be replaced with proper local index)
384 info.local_dune_to_local_hypre[local_ind] = 1;
385 } else {
386 // Mark as ghost/non-owned
387 info.local_dune_to_local_hypre[local_ind] = -1;
388 }
389 }
390
391 // STEP 2: Local Reordering & Owner-First Detection
392 // Create compact [0..N_owned-1] local HYPRE indexing for owned DOFs
393 // Simultaneously detect if owned DOFs appear before all ghost DOFs
394 bool owner_first = true;
395 bool visited_ghost = false; // Have we seen any ghost DOF yet?
396 std::size_t count = 0; // Counter for owned DOFs
397
398 for (std::size_t i = 0; i < info.local_dune_to_local_hypre.size(); ++i) {
399 if (info.local_dune_to_local_hypre[i] < 0) {
400 visited_ghost = true;
401 } else {
402 // This is an owned DOF - assign its local HYPRE index
403 info.local_dune_to_local_hypre[i] = count;
404 // Store the inverse mapping
405 info.local_hypre_to_local_dune.push_back(i);
406
407 // Check if we've seen ghost DOFs before this owner
408 owner_first = owner_first && !visited_ghost;
409 count += 1;
410 }
411 }
412
413 // Owner first need other copying of data
414 info.owner_first = owner_first;
415 info.N_owned = count;
416
417 // STEP 3: Global Offset Calculation
418 // Coordinate with other processes to determine global index ranges
419 // Each process owns a contiguous range of global indices
420 // Use HYPRE_Int to match the send buffer type exactly; on builds where
421 // HYPRE_Int is 64-bit (long long) using int here causes MPI_ERR_TRUNCATE.
422 std::vector<HYPRE_Int> dof_counts_per_process(collective_comm.size());
423 collective_comm.allgather(&info.N_owned, 1, dof_counts_per_process.data());
424
425 // Calculate this process's global offset (sum of DOFs in processes with lower rank)
426 info.dof_offset = std::accumulate(dof_counts_per_process.begin(),
427 dof_counts_per_process.begin() + collective_comm.rank(),
428 HYPRE_Int{0});
429
430 // STEP 4: Create Global Indices for Owned DOFs
431 // Convert local HYPRE indices to global HYPRE indices by adding offset
432 for (std::size_t i = 0; i < info.local_dune_to_local_hypre.size(); ++i) {
433 if (info.local_dune_to_local_hypre[i] >= 0) {
434 // Owned DOF: global index = local HYPRE index + this process's offset
436 } else {
437 info.local_dune_to_global_hypre[i] = -1;
438 }
439 }
440
441 if (collective_comm.rank() > 0) {
442 assert(info.dof_offset > 0);
443 }
444
445 // STEP 5: Exchange global indices for ghost DOFs
446 // After this call, ghost DOFs will have their correct global indices
447 comm.copyOwnerToAll(info.local_dune_to_global_hypre, info.local_dune_to_global_hypre);
448
449 return info;
450}
451
460template <typename MatrixType>
461SparsityPattern
462setupSparsityPattern(const MatrixType& matrix,
463 const ParallelInfo& par_info,
464 bool owner_first)
465{
466#if HYPRE_USING_CUDA || HYPRE_USING_HIP
468 return gpuistl::HypreInterface::setupSparsityPatternFromGpuMatrix(matrix, par_info, owner_first);
469 } else
470#endif
471 {
472 return setupSparsityPatternFromCpuMatrix(matrix, par_info, owner_first);
473 }
474}
475
484template <typename MatrixType>
485SparsityPattern
486setupSparsityPatternFromCpuMatrix(const MatrixType& matrix,
487 const ParallelInfo& par_info,
488 bool owner_first)
489{
490 SparsityPattern pattern;
491
492 // Determine the size for cols array based on owner_first
493 if (owner_first) {
494 std::size_t cols_size = 0;
495 // For owner_first=true case, we need to calculate how many owned entries there are
496 for (auto row = matrix.begin(); row != matrix.end(); ++row) {
497 const int rowIdx = row.index();
498 if (par_info.local_dune_to_local_hypre[rowIdx] >= 0) {
499 cols_size += row->size();
500 }
501 }
502 pattern.nnz = cols_size;
503 } else {
504 // Full matrix space case: all entries (including gaps)
505 pattern.nnz = matrix.nonzeroes();
506 }
507
508 // Setup host arrays
509 pattern.ncols.resize(par_info.N_owned);
510 pattern.rows.resize(par_info.N_owned);
511 pattern.cols.resize(pattern.nnz);
512
513 int pos = 0;
514 for (auto row = matrix.begin(); row != matrix.end(); ++row) {
515 const int rind = row.index();
516 const int local_rowIdx = par_info.local_dune_to_local_hypre[rind];
517
518 // For owner_first=true: skip ghost rows entirely
519 // For owner_first=false: process all rows (owned + ghost)
520 if (owner_first && local_rowIdx < 0) {
521 continue;
522 }
523
524 if (local_rowIdx >= 0) {
525 // This is an owned row - record its metadata
526 const int global_rowIdx = par_info.local_dune_to_global_hypre[rind];
527 pattern.rows[local_rowIdx] = global_rowIdx;
528 pattern.ncols[local_rowIdx] = row->size();
529 }
530
531 // Add column indices for this row
532 for (auto col = row->begin(); col != row->end(); ++col) {
533 const int global_colIdx = par_info.local_dune_to_global_hypre[col.index()];
534 assert(global_colIdx >= 0);
535 pattern.cols[pos++] = global_colIdx;
536 }
537 }
538
539 return pattern;
540}
541
555template <typename MatrixType>
556std::vector<HYPRE_Int>
557computeRowIndexes(const MatrixType& matrix,
558 const std::vector<HYPRE_Int>& ncols,
559 const std::vector<int>& local_dune_to_local_hypre,
560 bool owner_first)
561{
562 if (owner_first) {
563 // Simple contiguous case: prefix sum of ncols
564 std::vector<HYPRE_Int> row_indexes(ncols.size());
565 row_indexes[0] = 0;
566 for (std::size_t i = 1; i < ncols.size(); ++i) {
567 row_indexes[i] = row_indexes[i - 1] + ncols[i - 1];
568 }
569 return row_indexes;
570 } else {
571 // We need to compute the row indexes with mapping since we have gaps in the data
572#if HYPRE_USING_CUDA || HYPRE_USING_HIP
574 return gpuistl::HypreInterface::computeRowIndexesWithMappingGpu(matrix, local_dune_to_local_hypre);
575 } else
576#endif // HYPRE_USING_CUDA || HYPRE_USING_HIP
577 {
578 return computeRowIndexesWithMappingCpu(matrix, local_dune_to_local_hypre);
579 }
580 }
581}
582
602template <typename MatrixType>
603std::vector<HYPRE_Int>
604computeRowIndexesWithMappingCpu(const MatrixType& matrix, const std::vector<int>& local_dune_to_local_hypre)
605{
606 const int N = std::count_if(
607 local_dune_to_local_hypre.begin(), local_dune_to_local_hypre.end(), [](int val) { return val >= 0; });
608 std::vector<HYPRE_Int> row_indexes(N);
609 int data_position = 0; // Current position in the FULL (including ghost) matrix data
610 // Manually compute row starting positions by iterating through all matrix rows
611 for (auto row = matrix.begin(); row != matrix.end(); ++row) {
612 const int dune_row_idx = row.index();
613 const int hypre_row_idx = local_dune_to_local_hypre[dune_row_idx];
614
615 if (hypre_row_idx >= 0) {
616 // This is an owned row - record where its data starts in the FULL matrix
617 // Use hypre_row_idx as index (maps to Hypre ordering up to N_owned)
618 row_indexes[hypre_row_idx] = data_position;
619 }
620 // Always advance position (including ghost rows - this creates gaps for owned-only access)
621 data_position += row->size();
622 }
623 return row_indexes;
624}
625
626} // namespace Opm::gpuistl::HypreInterface
627
628#endif // OPM_HYPRE_SETUP_HPP
#define OPM_HYPRE_SAFE_CALL(expr)
Macro to wrap Hypre function calls with error checking.
Definition: HypreErrorHandling.hpp:96
#define HYPRE_SAFE_CALL(expr)
Short form macro for Hypre function calls (for backward compatibility)
Definition: HypreErrorHandling.hpp:102
Hierarchical collection of key/value pairs.
Definition: PropertyTree.hpp:39
T get(const std::string &key) const
std::vector< HYPRE_Int > computeRowIndexesWithMappingGpu(const GpuSparseMatrixWrapper< T, ForceLegacy > &gpu_matrix, const std::vector< int > &local_dune_to_local_hypre)
Compute row indexes for GPU matrix with ownership mapping.
Definition: gpuistl/hypreinterface/HypreSetup.hpp:150
linalg::HypreInterface::SparsityPattern setupSparsityPatternFromGpuMatrix(const GpuSparseMatrixWrapper< T, ForceLegacy > &gpu_matrix, const linalg::HypreInterface::ParallelInfo &par_info, bool owner_first)
Setup sparsity pattern from GPU matrix (GpuSparseMatrix)
Definition: gpuistl/hypreinterface/HypreSetup.hpp:71
Unified interface for Hypre operations with both CPU and GPU data structures.
Definition: hypreinterface/HypreCpuTransfers.hpp:36
SparsityPattern setupSparsityPatternFromCpuMatrix(const MatrixType &matrix, const ParallelInfo &par_info, bool owner_first)
Setup sparsity pattern from CPU matrix (BCRSMatrix)
Definition: hypreinterface/HypreSetup.hpp:486
void initialize(bool use_gpu_backend)
Initialize the Hypre library and set memory/execution policy.
Definition: hypreinterface/HypreSetup.hpp:84
HYPRE_IJVector createVector(HYPRE_Int N, HYPRE_Int dof_offset, const CommType &comm)
Create Hypre vector.
Definition: hypreinterface/HypreSetup.hpp:192
ParallelInfo setupHypreParallelInfoSerial(HYPRE_Int N)
Setup parallel information for Hypre in serial case.
Definition: hypreinterface/HypreSetup.hpp:274
ParallelInfo setupHypreParallelInfoParallel(const CommType &comm, const MatrixType &matrix)
Create mappings between Dune and HYPRE indexing for parallel decomposition.
Definition: hypreinterface/HypreSetup.hpp:358
void destroySolver(HYPRE_Solver solver)
Destroy Hypre solver.
Definition: hypreinterface/HypreSetup.hpp:214
SparsityPattern setupSparsityPattern(const MatrixType &matrix, const ParallelInfo &par_info, bool owner_first)
Setup sparsity pattern from matrix (automatically detects CPU/GPU type)
Definition: hypreinterface/HypreSetup.hpp:462
ParallelInfo setupHypreParallelInfo(const CommType &comm, const MatrixType &matrix)
Setup parallel information for Hypre (automatically detects serial/parallel)
Definition: hypreinterface/HypreSetup.hpp:258
void destroyMatrix(HYPRE_IJMatrix matrix)
Destroy Hypre matrix.
Definition: hypreinterface/HypreSetup.hpp:228
HYPRE_IJMatrix createMatrix(HYPRE_Int N, HYPRE_Int dof_offset, const CommType &comm)
Create Hypre matrix.
Definition: hypreinterface/HypreSetup.hpp:165
std::vector< HYPRE_Int > computeRowIndexesWithMappingCpu(const MatrixType &matrix, const std::vector< HYPRE_Int > &ncols, const std::vector< int > &local_dune_to_local_hypre, bool owner_first)
void setSolverParameters(HYPRE_Solver solver, const PropertyTree &prm, bool use_gpu_backend)
Set solver parameters from property tree.
Definition: hypreinterface/HypreSetup.hpp:128
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: hypreinterface/HypreSetup.hpp:557
void destroyVector(HYPRE_IJVector vector)
Destroy Hypre vector.
Definition: hypreinterface/HypreSetup.hpp:242
HYPRE_Solver createAMGSolver()
Create Hypre solver (BoomerAMG)
Definition: hypreinterface/HypreSetup.hpp:112
Type trait to detect if a type is a GPU type.
Definition: gpu_type_detection.hpp:40
Parallel domain decomposition information for HYPRE-Dune interface.
Definition: HypreDataStructures.hpp:37
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
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
HYPRE_Int N_owned
Number of DOFs owned by this MPI process.
Definition: HypreDataStructures.hpp:62
bool owner_first
Whether owned DOFs appear first in local Dune ordering.
Definition: HypreDataStructures.hpp:77
Compressed Sparse Row (CSR) sparsity pattern for HYPRE matrix assembly.
Definition: HypreDataStructures.hpp:86
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
std::vector< HYPRE_Int > ncols
Non-zero entries per owned row (size: N_owned)
Definition: HypreDataStructures.hpp:88
std::vector< HYPRE_BigInt > rows
Global row indices for owned rows (size: N_owned)
Definition: HypreDataStructures.hpp:91