Opm::gpuistl::HypreInterface Namespace Reference

Unified interface for Hypre operations with both CPU and GPU data structures. More...

Typedefs

using HostArrays = ::Opm::gpuistl::HypreHostDataArrays
 
using DeviceArrays = ::Opm::gpuistl::HypreDeviceDataArrays
 
using ParallelInfo = ::Opm::gpuistl::ParallelInfo
 
using SparsityPattern = ::Opm::gpuistl::SparsityPattern
 

Functions

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 >
SparsityPattern setupSparsityPatternFromGpuMatrix (const GpuSparseMatrix< T > &gpu_matrix, const ParallelInfo &par_info, bool owner_first)
 
template<typename T >
std::vector< HYPRE_Int > computeRowIndexesWithMappingGpu (const GpuSparseMatrix< T > &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...
 

Detailed Description

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:

  1. Input type is CPU and backend acceleration is CPU
  2. Input type is CPU and backend acceleration is GPU
  3. Input type is GPU and backend acceleration is GPU
  4. 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.

Typedef Documentation

◆ DeviceArrays

◆ HostArrays

◆ ParallelInfo

◆ SparsityPattern

Function Documentation

◆ 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
matrixThe matrix to analyze
ncolsArray containing number of columns per row
local_dune_to_local_hypreMapping from Dune indices to Hypre indices
owner_firstWhether 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 
)

Referenced by computeRowIndexes().

◆ 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
matrixMatrix to analyze
local_dune_to_local_hypreMapping 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 >
std::vector< HYPRE_Int > Opm::gpuistl::HypreInterface::computeRowIndexesWithMappingGpu ( const GpuSparseMatrix< T > &  gpu_matrix,
const std::vector< int > &  local_dune_to_local_hypre 
)

Referenced by computeRowIndexes().

◆ createAMGSolver()

HYPRE_Solver Opm::gpuistl::HypreInterface::createAMGSolver ( )
inline

Create Hypre solver (BoomerAMG)

Create Hypre BoomerAMG solver.

Returns
HYPRE_Solver The created solver handle
Exceptions
HypreErrorif solver creation fails

References OPM_HYPRE_SAFE_CALL.

Referenced by Hypre::HyprePreconditioner< M, X, Y, Comm >::HyprePreconditioner().

◆ createMatrix()

template<typename CommType >
HYPRE_IJMatrix Opm::gpuistl::HypreInterface::createMatrix ( HYPRE_Int  N,
HYPRE_Int  dof_offset,
const CommType &  comm 
)

Create Hypre matrix.

Parameters
NNumber of rows/columns
dof_offsetOffset for the matrix
commThe communicator
Returns
HYPRE_IJMatrix The created matrix handle
Exceptions
HypreErrorif matrix creation fails

References OPM_HYPRE_SAFE_CALL.

Referenced by Opm::getMatrixRowColoring(), and Hypre::HyprePreconditioner< M, X, Y, Comm >::HyprePreconditioner().

◆ createVector()

template<typename CommType >
HYPRE_IJVector Opm::gpuistl::HypreInterface::createVector ( HYPRE_Int  N,
HYPRE_Int  dof_offset,
const CommType &  comm 
)

Create Hypre vector.

Parameters
NSize of the vector
dof_offsetOffset for the vector
commThe communicator
Returns
HYPRE_IJVector The created vector handle
Exceptions
HypreErrorif vector creation fails

References OPM_HYPRE_SAFE_CALL.

Referenced by Hypre::HyprePreconditioner< M, X, Y, Comm >::HyprePreconditioner().

◆ destroyMatrix()

void Opm::gpuistl::HypreInterface::destroyMatrix ( HYPRE_IJMatrix  matrix)
inline

Destroy Hypre matrix.

Parameters
matrixThe matrix handle to destroy
Exceptions
HypreErrorif matrix destruction fails

References OPM_HYPRE_SAFE_CALL.

Referenced by Hypre::HyprePreconditioner< M, X, Y, Comm >::~HyprePreconditioner().

◆ destroySolver()

void Opm::gpuistl::HypreInterface::destroySolver ( HYPRE_Solver  solver)
inline

Destroy Hypre solver.

Parameters
solverThe solver handle to destroy
Exceptions
HypreErrorif solver destruction fails

References OPM_HYPRE_SAFE_CALL.

Referenced by Hypre::HyprePreconditioner< M, X, Y, Comm >::~HyprePreconditioner().

◆ destroyVector()

void Opm::gpuistl::HypreInterface::destroyVector ( HYPRE_IJVector  vector)
inline

Destroy Hypre vector.

Parameters
vectorThe vector handle to destroy
Exceptions
HypreErrorif vector destruction fails

References OPM_HYPRE_SAFE_CALL.

Referenced by Hypre::HyprePreconditioner< M, X, Y, Comm >::~HyprePreconditioner().

◆ 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_matrixThe Hypre matrix to read values from
ncolsArray containing number of columns per row
rowsArray containing row indices
colsArray containing column indices
use_gpu_backendWhether the Hypre backend is using GPU acceleration
Returns
Vector containing the retrieved matrix values
Exceptions
HypreErrorif the operation fails
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

◆ 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 
)

Extract owned vector values in the order expected by HYPRE.

Referenced by transferCpuVectorToHypre().

◆ 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 
)

Distribute HYPRE vector values back to original vector positions.

Referenced by transferHypreToCpuVector().

◆ 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

Set solver parameters from property tree.

Parameters
solverThe solver handle
prmProperty tree containing configuration parameters
use_gpu_backendWhether using GPU backend
Exceptions
HypreErrorif parameter setting fails

References Opm::PropertyTree::get(), and HYPRE_SAFE_CALL.

Referenced by Hypre::HyprePreconditioner< M, X, Y, Comm >::HyprePreconditioner().

◆ setupHypreParallelInfo()

template<typename CommType , typename MatrixType >
ParallelInfo Opm::gpuistl::HypreInterface::setupHypreParallelInfo ( const CommType &  comm,
const MatrixType &  matrix 
)

Setup parallel information for Hypre (automatically detects serial/parallel)

Parameters
commThe communication object (serial or parallel)
matrixThe matrix to analyze
Returns
ParallelInfo structure containing mappings and offsets

References setupHypreParallelInfoParallel(), and setupHypreParallelInfoSerial().

Referenced by Hypre::HyprePreconditioner< M, X, Y, Comm >::HyprePreconditioner().

◆ 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:

  1. 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
  2. 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
  3. 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:

  1. Ownership Detection: Scan Dune's index set to identify owner vs ghost DOFs
  2. Local Reordering: Create compact [0..N_owned-1] indexing for owned DOFs
  3. Owner-First Detection: Determine if all owned DOFs appear before ghost DOFs
  4. Global Offset Calculation: Coordinate with other processes to assign global ranges
  5. 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
commDune communication object with parallel index set
matrixDistributed 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 
)

Setup sparsity pattern from matrix (automatically detects CPU/GPU type)

Parameters
matrixThe matrix to analyze (CPU or GPU type)
par_infoParallel information structure
owner_firstWhether all owned DOFs come first in the ordering
Returns
Sparsity pattern information

References setupSparsityPatternFromCpuMatrix(), and setupSparsityPatternFromGpuMatrix().

Referenced by Hypre::HyprePreconditioner< M, X, Y, Comm >::HyprePreconditioner().

◆ setupSparsityPatternFromCpuMatrix()

template<typename MatrixType >
SparsityPattern Opm::gpuistl::HypreInterface::setupSparsityPatternFromCpuMatrix ( const MatrixType &  matrix,
const ParallelInfo par_info,
bool  owner_first 
)

Setup sparsity pattern from CPU matrix (BCRSMatrix)

Parameters
matrixThe matrix to analyze
par_infoParallel information structure
owner_firstWhether all owned DOFs come first in the ordering
Returns
Sparsity pattern information

References Opm::gpuistl::SparsityPattern::cols, Opm::gpuistl::ParallelInfo::local_dune_to_global_hypre, Opm::gpuistl::ParallelInfo::local_dune_to_local_hypre, Opm::gpuistl::ParallelInfo::N_owned, Opm::gpuistl::SparsityPattern::ncols, Opm::gpuistl::SparsityPattern::nnz, and Opm::gpuistl::SparsityPattern::rows.

Referenced by setupSparsityPattern().

◆ setupSparsityPatternFromGpuMatrix()

template<typename T >
SparsityPattern Opm::gpuistl::HypreInterface::setupSparsityPatternFromGpuMatrix ( const GpuSparseMatrix< T > &  gpu_matrix,
const ParallelInfo par_info,
bool  owner_first 
)

Referenced by setupSparsityPattern().

◆ transferCpuVectorToHypre()

template<typename VectorType >
void Opm::gpuistl::HypreInterface::transferCpuVectorToHypre ( const VectorType &  cpu_vec,
HYPRE_IJVector  hypre_vec,
HypreHostDataArrays host_arrays,
const HypreDeviceDataArrays device_arrays,
const ParallelInfo par_info,
bool  use_gpu_backend 
)

◆ transferHypreToCpuVector()

template<typename VectorType >
void Opm::gpuistl::HypreInterface::transferHypreToCpuVector ( HYPRE_IJVector  hypre_vec,
VectorType &  cpu_vec,
HypreHostDataArrays host_arrays,
const HypreDeviceDataArrays device_arrays,
const ParallelInfo par_info,
bool  use_gpu_backend 
)

◆ 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 
)

Transfer vector from Hypre to any vector type (CPU or GPU)

References transferHypreToCpuVector().

Referenced by Hypre::HyprePreconditioner< M, X, Y, Comm >::apply().

◆ 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 
)

Transfer vector to Hypre from any vector type (CPU or GPU)

References transferCpuVectorToHypre().

Referenced by Hypre::HyprePreconditioner< M, X, Y, Comm >::apply().

◆ updateMatrixFromCpuMatrix()

template<typename MatrixType >
void Opm::gpuistl::HypreInterface::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 
)

◆ 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 
)