Unified interface for Hypre operations with both CPU and GPU data structures.
More...
|
| void | initialize (bool use_gpu_backend) |
| | Initialize the Hypre library and set memory/execution policy. More...
|
| |
| HYPRE_Solver | createAMGSolver () |
| | Create Hypre solver (BoomerAMG) More...
|
| |
| void | setSolverParameters (HYPRE_Solver solver, const PropertyTree &prm, bool use_gpu_backend) |
| | Set solver parameters from property tree. More...
|
| |
| template<typename CommType > |
| HYPRE_IJMatrix | createMatrix (HYPRE_Int N, HYPRE_Int dof_offset, const CommType &comm) |
| | Create Hypre matrix. More...
|
| |
| template<typename CommType > |
| HYPRE_IJVector | createVector (HYPRE_Int N, HYPRE_Int dof_offset, const CommType &comm) |
| | Create Hypre vector. More...
|
| |
| void | destroySolver (HYPRE_Solver solver) |
| | Destroy Hypre solver. More...
|
| |
| void | destroyMatrix (HYPRE_IJMatrix matrix) |
| | Destroy Hypre matrix. More...
|
| |
| void | destroyVector (HYPRE_IJVector vector) |
| | Destroy Hypre vector. More...
|
| |
| template<typename CommType , typename MatrixType > |
| ParallelInfo | setupHypreParallelInfo (const CommType &comm, const MatrixType &matrix) |
| | Setup parallel information for Hypre (automatically detects serial/parallel) More...
|
| |
| template<typename MatrixType > |
| SparsityPattern | setupSparsityPattern (const MatrixType &matrix, const ParallelInfo &par_info, bool owner_first) |
| | Setup sparsity pattern from matrix (automatically detects CPU/GPU type) More...
|
| |
| template<typename MatrixType > |
| 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. More...
|
| |
| template<typename VectorType > |
| void | transferVectorToHypre (const VectorType &vec, HYPRE_IJVector hypre_vec, HostArrays &host_arrays, const DeviceArrays &device_arrays, const ParallelInfo &par_info, bool use_gpu_backend) |
| | Transfer vector to Hypre from any vector type (CPU or GPU) More...
|
| |
| template<typename VectorType > |
| void | transferVectorFromHypre (HYPRE_IJVector hypre_vec, VectorType &vec, HostArrays &host_arrays, const DeviceArrays &device_arrays, const ParallelInfo &par_info, bool use_gpu_backend) |
| | Transfer vector from Hypre to any vector type (CPU or GPU) More...
|
| |
| template<typename MatrixType > |
| void | updateMatrixValues (const MatrixType &matrix, HYPRE_IJMatrix hypre_matrix, const SparsityPattern &sparsity_pattern, const HostArrays &host_arrays, const DeviceArrays &device_arrays, bool use_gpu_backend) |
| | Update matrix values in Hypre. More...
|
| |
| template<typename VectorType > |
| void | setContinuousVectorForHypre (const VectorType &v, std::vector< HYPRE_Real > &continuous_vector_values, const std::vector< int > &local_hypre_to_local_dune) |
| | Extract owned vector values in the order expected by HYPRE. More...
|
| |
| template<typename VectorType > |
| void | setDuneVectorFromContinuousVector (VectorType &v, const std::vector< HYPRE_Real > &continuous_vector_values, const std::vector< int > &local_hypre_to_local_dune) |
| | Distribute HYPRE vector values back to original vector positions. More...
|
| |
| template<typename VectorType > |
| void | transferCpuVectorToHypre (const VectorType &cpu_vec, HYPRE_IJVector hypre_vec, HypreHostDataArrays &host_arrays, const HypreDeviceDataArrays &device_arrays, const ParallelInfo &par_info, bool use_gpu_backend) |
| | Transfer CPU vector to Hypre vector. More...
|
| |
| template<typename VectorType > |
| void | transferHypreToCpuVector (HYPRE_IJVector hypre_vec, VectorType &cpu_vec, HypreHostDataArrays &host_arrays, const HypreDeviceDataArrays &device_arrays, const ParallelInfo &par_info, bool use_gpu_backend) |
| | Transfer Hypre vector to CPU vector. More...
|
| |
| template<typename MatrixType > |
| void | updateMatrixFromCpuMatrix (const MatrixType &cpu_matrix, HYPRE_IJMatrix hypre_matrix, const SparsityPattern &sparsity_pattern, const HypreHostDataArrays &host_arrays, const HypreDeviceDataArrays &device_arrays, bool use_gpu_backend) |
| | Update Hypre matrix from CPU matrix Uses HYPRE_IJMatrixSetValues2 with pre-computed row_indexes, which allows us to use the original CPU matrix data (with potential ghost values) directly. More...
|
| |
| template<typename VectorType > |
| void | setContinuousGpuVectorForHypre (const VectorType &v, std::vector< HYPRE_Real > &continuous_vector_values, const std::vector< int > &local_hypre_to_local_dune) |
| |
| template<typename VectorType > |
| void | setGpuVectorFromContinuousVector (VectorType &v, const std::vector< HYPRE_Real > &continuous_vector_values, const std::vector< int > &local_hypre_to_local_dune) |
| |
| template<typename T , bool ForceLegacy> |
| SparsityPattern | setupSparsityPatternFromGpuMatrix (const GpuSparseMatrixWrapper< T, ForceLegacy > &gpu_matrix, const ParallelInfo &par_info, bool owner_first) |
| |
| template<typename T , bool ForceLegacy> |
| std::vector< HYPRE_Int > | computeRowIndexesWithMappingGpu (const GpuSparseMatrixWrapper< T, ForceLegacy > &gpu_matrix, const std::vector< int > &local_dune_to_local_hypre) |
| |
| ParallelInfo | setupHypreParallelInfoSerial (HYPRE_Int N) |
| | Setup parallel information for Hypre in serial case. More...
|
| |
| template<typename CommType , typename MatrixType > |
| ParallelInfo | setupHypreParallelInfoParallel (const CommType &comm, const MatrixType &matrix) |
| | Create mappings between Dune and HYPRE indexing for parallel decomposition. More...
|
| |
| template<typename MatrixType > |
| SparsityPattern | setupSparsityPatternFromCpuMatrix (const MatrixType &matrix, const ParallelInfo &par_info, bool owner_first) |
| | Setup sparsity pattern from CPU matrix (BCRSMatrix) More...
|
| |
| template<typename MatrixType > |
| 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) |
| |
| template<typename MatrixType > |
| std::vector< HYPRE_Int > | computeRowIndexesWithMappingCpu (const MatrixType &matrix, const std::vector< int > &local_dune_to_local_hypre) |
| | Compute row indexes for CPU matrix with ownership mapping. More...
|
| |
| std::vector< HYPRE_Real > | getMatrixValues (HYPRE_IJMatrix hypre_matrix, const std::vector< HYPRE_Int > &ncols, const std::vector< HYPRE_BigInt > &rows, const std::vector< HYPRE_BigInt > &cols, bool use_gpu_backend=false) |
| | Get matrix values from Hypre matrix. More...
|
| |
Unified interface for Hypre operations with both CPU and GPU data structures.
This namespace provides utilities for working with Hypre resources and transferring data between CPU/GPU data structures and Hypre handles. It handles type detection and automatically chooses the most efficient transfer method:
- For GPU types: Zero-copy device-to-device transfers when using GPU backend
- For CPU types: Host-to-device transfers when using GPU backend
- For CPU types with CPU backend: Direct host memory usage
Supports four use cases:
- Input type is CPU and backend acceleration is CPU
- Input type is CPU and backend acceleration is GPU
- Input type is GPU and backend acceleration is GPU
- Input type is GPU and backend acceleration is CPU
- Note
- Error handling: All functions throw HypreError exceptions when Hypre operations fail.
-
This is a consolidated version that includes all functionality across multiple files.
◆ DeviceArrays
◆ HostArrays
◆ ParallelInfo
◆ SparsityPattern
◆ computeRowIndexes()
template<typename MatrixType >
| std::vector< HYPRE_Int > Opm::gpuistl::HypreInterface::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.
Chooses the appropriate row index computation method based on owner_first and matrix type.
- For owner_first=true: Simple prefix sum for contiguous data
- For owner_first=false: Mapping-based approach to handle gaps from ghost rows
- Parameters
-
| matrix | The matrix to analyze |
| ncols | Array containing number of columns per row |
| local_dune_to_local_hypre | Mapping from Dune indices to Hypre indices |
| owner_first | Whether all owned DOFs come first in the ordering |
- Returns
- Vector containing row indexes for HYPRE_IJMatrixSetValues2
References computeRowIndexesWithMappingCpu(), and computeRowIndexesWithMappingGpu().
Referenced by Hypre::HyprePreconditioner< M, X, Y, Comm >::HyprePreconditioner().
◆ computeRowIndexesWithMappingCpu() [1/2]
template<typename MatrixType >
| std::vector< HYPRE_Int > Opm::gpuistl::HypreInterface::computeRowIndexesWithMappingCpu |
( |
const MatrixType & |
matrix, |
|
|
const std::vector< HYPRE_Int > & |
ncols, |
|
|
const std::vector< int > & |
local_dune_to_local_hypre, |
|
|
bool |
owner_first |
|
) |
| |
◆ computeRowIndexesWithMappingCpu() [2/2]
template<typename MatrixType >
| std::vector< HYPRE_Int > Opm::gpuistl::HypreInterface::computeRowIndexesWithMappingCpu |
( |
const MatrixType & |
matrix, |
|
|
const std::vector< int > & |
local_dune_to_local_hypre |
|
) |
| |
Compute row indexes for CPU matrix with ownership mapping.
Creates row_indexes that point directly into the FULL matrix data, enabling using the original matrix data without any data copying. This is necessary for some cases where ghost rows can come before all owned rows.
Example with local_dune_to_local_hypre = [-1, 0, -1, 1]: Original matrix data (with gaps for ghost rows): Row 0 (ghost): {_, _} at positions 0-1 (skipped) Row 1 (owned): {1.0, 2.0} at positions 2-3 Row 2 (ghost): {_} at position 4 (skipped) Row 3 (owned): {3.0} at position 5 Resulting row_indexes: {2, 5} This tells HYPRE that Row 0 data starts at position 2, Row 1 data starts at position 5
- Parameters
-
| matrix | Matrix to analyze |
| local_dune_to_local_hypre | Mapping from Dune indices to Hypre indices (-1 for ghost) |
- Returns
- Vector containing row indexes that map to full matrix data positions
◆ computeRowIndexesWithMappingGpu()
template<typename T , bool ForceLegacy>
| std::vector< HYPRE_Int > Opm::gpuistl::HypreInterface::computeRowIndexesWithMappingGpu |
( |
const GpuSparseMatrixWrapper< T, ForceLegacy > & |
gpu_matrix, |
|
|
const std::vector< int > & |
local_dune_to_local_hypre |
|
) |
| |
◆ createAMGSolver()
| HYPRE_Solver Opm::gpuistl::HypreInterface::createAMGSolver |
( |
| ) |
|
|
inline |
◆ createMatrix()
template<typename CommType >
| HYPRE_IJMatrix Opm::gpuistl::HypreInterface::createMatrix |
( |
HYPRE_Int |
N, |
|
|
HYPRE_Int |
dof_offset, |
|
|
const CommType & |
comm |
|
) |
| |
◆ createVector()
template<typename CommType >
| HYPRE_IJVector Opm::gpuistl::HypreInterface::createVector |
( |
HYPRE_Int |
N, |
|
|
HYPRE_Int |
dof_offset, |
|
|
const CommType & |
comm |
|
) |
| |
◆ destroyMatrix()
| void Opm::gpuistl::HypreInterface::destroyMatrix |
( |
HYPRE_IJMatrix |
matrix | ) |
|
|
inline |
◆ destroySolver()
| void Opm::gpuistl::HypreInterface::destroySolver |
( |
HYPRE_Solver |
solver | ) |
|
|
inline |
◆ destroyVector()
| void Opm::gpuistl::HypreInterface::destroyVector |
( |
HYPRE_IJVector |
vector | ) |
|
|
inline |
◆ getMatrixValues()
| std::vector< HYPRE_Real > Opm::gpuistl::HypreInterface::getMatrixValues |
( |
HYPRE_IJMatrix |
hypre_matrix, |
|
|
const std::vector< HYPRE_Int > & |
ncols, |
|
|
const std::vector< HYPRE_BigInt > & |
rows, |
|
|
const std::vector< HYPRE_BigInt > & |
cols, |
|
|
bool |
use_gpu_backend = false |
|
) |
| |
|
inline |
Get matrix values from Hypre matrix.
Retrieves the values from a Hypre matrix for verification or debugging purposes.
- Parameters
-
| hypre_matrix | The Hypre matrix to read values from |
| ncols | Array containing number of columns per row |
| rows | Array containing row indices |
| cols | Array containing column indices |
| use_gpu_backend | Whether the Hypre backend is using GPU acceleration |
- Returns
- Vector containing the retrieved matrix values
- Exceptions
-
- Note
- For GPU backend, uses HYPRE_ParCSRMatrixGetRow since HYPRE_IJMatrixGetValues not supported on GPU
References OPM_HYPRE_SAFE_CALL.
◆ initialize()
| void Opm::gpuistl::HypreInterface::initialize |
( |
bool |
use_gpu_backend | ) |
|
|
inline |
Initialize the Hypre library and set memory/execution policy.
- Parameters
-
| use_gpu_backend | Whether to use GPU backend acceleration |
- Exceptions
-
References OPM_HYPRE_SAFE_CALL.
Referenced by Dune::CartesianIndexMapper< Dune::ALUGrid< 3, 3, Dune::cube, Dune::nonconforming > >::GlobalIndexDataHandle::GlobalIndexDataHandle(), Hypre::HyprePreconditioner< M, X, Y, Comm >::HyprePreconditioner(), Opm::CopyRestrictProlong< Grid, Container >::prolongLocal(), and Opm::CopyRestrictProlong< Grid, Container >::restrictLocal().
◆ setContinuousGpuVectorForHypre()
template<typename VectorType >
| void Opm::gpuistl::HypreInterface::setContinuousGpuVectorForHypre |
( |
const VectorType & |
v, |
|
|
std::vector< HYPRE_Real > & |
continuous_vector_values, |
|
|
const std::vector< int > & |
local_hypre_to_local_dune |
|
) |
| |
◆ setContinuousVectorForHypre()
template<typename VectorType >
| void Opm::gpuistl::HypreInterface::setContinuousVectorForHypre |
( |
const VectorType & |
v, |
|
|
std::vector< HYPRE_Real > & |
continuous_vector_values, |
|
|
const std::vector< int > & |
local_hypre_to_local_dune |
|
) |
| |
◆ setDuneVectorFromContinuousVector()
template<typename VectorType >
| void Opm::gpuistl::HypreInterface::setDuneVectorFromContinuousVector |
( |
VectorType & |
v, |
|
|
const std::vector< HYPRE_Real > & |
continuous_vector_values, |
|
|
const std::vector< int > & |
local_hypre_to_local_dune |
|
) |
| |
◆ setGpuVectorFromContinuousVector()
template<typename VectorType >
| void Opm::gpuistl::HypreInterface::setGpuVectorFromContinuousVector |
( |
VectorType & |
v, |
|
|
const std::vector< HYPRE_Real > & |
continuous_vector_values, |
|
|
const std::vector< int > & |
local_hypre_to_local_dune |
|
) |
| |
◆ setSolverParameters()
| void Opm::gpuistl::HypreInterface::setSolverParameters |
( |
HYPRE_Solver |
solver, |
|
|
const PropertyTree & |
prm, |
|
|
bool |
use_gpu_backend |
|
) |
| |
|
inline |
◆ setupHypreParallelInfo()
template<typename CommType , typename MatrixType >
| ParallelInfo Opm::gpuistl::HypreInterface::setupHypreParallelInfo |
( |
const CommType & |
comm, |
|
|
const MatrixType & |
matrix |
|
) |
| |
◆ setupHypreParallelInfoParallel()
template<typename CommType , typename MatrixType >
| ParallelInfo Opm::gpuistl::HypreInterface::setupHypreParallelInfoParallel |
( |
const CommType & |
comm, |
|
|
const MatrixType & |
matrix |
|
) |
| |
|
inline |
Create mappings between Dune and HYPRE indexing for parallel decomposition.
This function interfaces between Dune's distributed matrix representation and HYPRE's global indexing requirements. Note that Dune uses local indices with owner/ghost classification, while HYPRE requires globally consistent indices for all DOFs.
Index System Overview:
Dune Local Indices: Each MPI process has local DOFs indexed 0..N_local-1, where:
- "Owner" DOFs: This process is responsible for their values/updates
- "Ghost" DOFs: Copies of DOFs owned by other processes (for stencil access)
- These can be interleaved or arranged in "owner-first" order where all owned DOFs come first
HYPRE Global Indices: All DOFs across all processes have unique global indices:
- Process 0 owns global indices [0, N_owned_0-1]
- Process 1 owns global indices [N_owned_0, N_owned_0 + N_owned_1-1]
- Process k owns global indices [offset_k, offset_k + N_owned_k-1]
Index Mappings Created:
- local_dune_to_local_hypre: Maps Dune local index → HYPRE local index
- Size: N_local (all DOFs)
- Value: [0..N_owned-1] for owned DOFs, -1 for ghost DOFs
- Purpose: Identify owned DOFs and their compact local ordering
- local_dune_to_global_hypre: Maps Dune local index → HYPRE global index
- Size: N_local (all DOFs)
- Value: [offset..offset+N_owned-1] for owned, actual global index for ghost
- Purpose: Matrix assembly and vector operations with global indexing
- local_hypre_to_local_dune: Maps HYPRE local index → Dune local index
- Size: N_owned (owned DOFs only)
- Purpose: Transfer data from HYPRE back to Dune structures
Algorithm Steps:
- Ownership Detection: Scan Dune's index set to identify owner vs ghost DOFs
- Local Reordering: Create compact [0..N_owned-1] indexing for owned DOFs
- Owner-First Detection: Determine if all owned DOFs appear before ghost DOFs
- Global Offset Calculation: Coordinate with other processes to assign global ranges
- Ghost Communication: Exchange global indices for ghost DOFs
Example (2 processes, 3 DOFs each):
Process 0: Dune local [0,1,2] → owned=[0,1], ghost=[2 from P1]
- local_dune_to_local_hypre = [0, 1, -1]
- local_dune_to_global_hypre = [0, 1, 2] (after communication)
- local_hypre_to_local_dune = [0, 1]
- N_owned=2, dof_offset=0, owner_first=true
Process 1: Dune local [0,1,2] → owned=[1,2], ghost=[0 from P0]
- local_dune_to_local_hypre = [-1, 0, 1]
- local_dune_to_global_hypre = [0, 2, 3] (after communication)
- local_hypre_to_local_dune = [1, 2]
- N_owned=2, dof_offset=2, owner_first=false
- Parameters
-
| comm | Dune communication object with parallel index set |
| matrix | Distributed matrix for consistency checks |
- Returns
- ParallelInfo with all mappings, offsets, and metadata
- Note
- May modify comm if index set holes are detected (rare in OPM)
References Opm::gpuistl::ParallelInfo::dof_offset, Opm::gpuistl::ParallelInfo::local_dune_to_global_hypre, Opm::gpuistl::ParallelInfo::local_dune_to_local_hypre, Opm::gpuistl::ParallelInfo::local_hypre_to_local_dune, Opm::gpuistl::ParallelInfo::N_owned, and Opm::gpuistl::ParallelInfo::owner_first.
Referenced by setupHypreParallelInfo().
◆ setupHypreParallelInfoSerial()
| ParallelInfo Opm::gpuistl::HypreInterface::setupHypreParallelInfoSerial |
( |
HYPRE_Int |
N | ) |
|
|
inline |
◆ setupSparsityPattern()
template<typename MatrixType >
| SparsityPattern Opm::gpuistl::HypreInterface::setupSparsityPattern |
( |
const MatrixType & |
matrix, |
|
|
const ParallelInfo & |
par_info, |
|
|
bool |
owner_first |
|
) |
| |
◆ setupSparsityPatternFromCpuMatrix()
template<typename MatrixType >
| SparsityPattern Opm::gpuistl::HypreInterface::setupSparsityPatternFromCpuMatrix |
( |
const MatrixType & |
matrix, |
|
|
const ParallelInfo & |
par_info, |
|
|
bool |
owner_first |
|
) |
| |
◆ setupSparsityPatternFromGpuMatrix()
template<typename T , bool ForceLegacy>
◆ transferCpuVectorToHypre()
template<typename VectorType >
◆ transferHypreToCpuVector()
template<typename VectorType >
◆ transferVectorFromHypre()
template<typename VectorType >
| void Opm::gpuistl::HypreInterface::transferVectorFromHypre |
( |
HYPRE_IJVector |
hypre_vec, |
|
|
VectorType & |
vec, |
|
|
HostArrays & |
host_arrays, |
|
|
const DeviceArrays & |
device_arrays, |
|
|
const ParallelInfo & |
par_info, |
|
|
bool |
use_gpu_backend |
|
) |
| |
◆ transferVectorToHypre()
template<typename VectorType >
| void Opm::gpuistl::HypreInterface::transferVectorToHypre |
( |
const VectorType & |
vec, |
|
|
HYPRE_IJVector |
hypre_vec, |
|
|
HostArrays & |
host_arrays, |
|
|
const DeviceArrays & |
device_arrays, |
|
|
const ParallelInfo & |
par_info, |
|
|
bool |
use_gpu_backend |
|
) |
| |
◆ updateMatrixFromCpuMatrix()
template<typename MatrixType >
◆ updateMatrixValues()
template<typename MatrixType >
| void Opm::gpuistl::HypreInterface::updateMatrixValues |
( |
const MatrixType & |
matrix, |
|
|
HYPRE_IJMatrix |
hypre_matrix, |
|
|
const SparsityPattern & |
sparsity_pattern, |
|
|
const HostArrays & |
host_arrays, |
|
|
const DeviceArrays & |
device_arrays, |
|
|
bool |
use_gpu_backend |
|
) |
| |
|