Block-Structured AMR Software Framework
amrex::HypreSolver< MSS > Class Template Reference

Solve Ax = b using HYPRE's generic IJ matrix format where A is a sparse matrix specified using the compressed sparse row (CSR) format. More...

#include <AMReX_HypreSolver.H>

Public Member Functions

template<class Marker , class Filler >
 HypreSolver (Vector< IndexType > const &a_index_type, IntVect const &a_nghost, Geometry const &a_geom, BoxArray const &a_grids, DistributionMapping const &a_dmap, Marker &&a_marker, Filler &&a_filler, int a_verbose=0, std::string a_options_namespace="hypre")
 
template<class MF , std::enable_if_t< IsFabArray< MF >::value &&std::is_same_v< typename MF::value_type, HYPRE_Real >, int > = 0>
void solve (Vector< MF * > const &a_soln, Vector< MF const * > const &a_rhs, HYPRE_Real rel_tol, HYPRE_Real abs_tol, int max_iter)
 
int getNumIters () const
 
HYPRE_Real getFinalResidualNorm () const
 
HYPRE_IJMatrix getA () const
 
HYPRE_IJVector getb () const
 
HYPRE_IJVector getx () const
 
template<class Marker >
std::enable_if_t< IsCallable< Marker, int, int, int, int, int >::value > fill_local_id (Marker const &marker)
 
void fill_global_id ()
 
template<class Filler , std::enable_if_t< IsCallable< Filler, int, int, int, int, int, Array4< HYPRE_Int const > const *, HYPRE_Int &, HYPRE_Int *, HYPRE_Real * >::value, int > FOO = 0>
void fill_matrix (Filler const &filler)
 
template<class MF , std::enable_if_t< IsFabArray< MF >::value &&std::is_same_v< typename MF::value_type, HYPRE_Real >, int > = 0>
void load_vectors (Vector< MF * > const &a_soln, Vector< MF const * > const &a_rhs)
 
template<class MF , std::enable_if_t< IsFabArray< MF >::value &&std::is_same_v< typename MF::value_type, HYPRE_Real >, int > = 0>
void get_solution (Vector< MF * > const &a_soln)
 

Private Attributes

int m_nvars
 
Vector< IndexTypem_index_type
 
IntVect m_nghost
 
Geometry m_geom
 
Vector< BoxArraym_grids
 
DistributionMapping m_dmap
 
int m_verbose
 
std::string m_options_namespace
 
MPI_Comm m_comm = MPI_COMM_NULL
 
Vector< std::unique_ptr< iMultiFab > > m_owner_mask
 
Vector< iMultiFabm_local_id
 
Vector< FabArray< BaseFab< HYPRE_Int > > > m_global_id
 
LayoutData< Gpu::DeviceVector< HYPRE_Int > > m_global_id_vec
 
LayoutData< Gpu::DeviceVector< int > > m_cell_offset
 
Vector< LayoutData< HYPRE_Int > > m_nrows_grid
 
Vector< LayoutData< HYPRE_Int > > m_id_offset
 
LayoutData< HYPRE_Int > m_nrows
 
HYPRE_Int m_nrows_proc
 
std::unique_ptr< HypreIJIfacem_hypre_ij
 
HYPRE_IJMatrix m_A = nullptr
 
HYPRE_IJVector m_b = nullptr
 
HYPRE_IJVector m_x = nullptr
 

Detailed Description

template<int MSS>
class amrex::HypreSolver< MSS >

Solve Ax = b using HYPRE's generic IJ matrix format where A is a sparse matrix specified using the compressed sparse row (CSR) format.

An example of using HypreSolver to solve Poisson's equation is located in Tests/LinearSolvers/Hypre

Constructor & Destructor Documentation

◆ HypreSolver()

template<int MSS>
template<class Marker , class Filler >
amrex::HypreSolver< MSS >::HypreSolver ( Vector< IndexType > const &  a_index_type,
IntVect const &  a_nghost,
Geometry const &  a_geom,
BoxArray const &  a_grids,
DistributionMapping const &  a_dmap,
Marker &&  a_marker,
Filler &&  a_filler,
int  a_verbose = 0,
std::string  a_options_namespace = "hypre" 
)

Constructor

Parameters
[in]a_index_typenodality of the data
[in]a_nghostnumber of ghosts cells
[in]a_geomGeometry of problem domain
[in]a_gridsBoxArray of computational grids
[in]a_dmapDistributionMapping
[in]a_markerfunctor that returns whether the variable n at (i,j,k) in Box boxno (local index) is valid (i.e., not exactly on Dirichlet boundary). function signature:
(int boxno, int i, int j, int k, int n) -> bool
[in]a_fillerfunctor that fills the row in the matrix A for variable n at (i,j,k) in Box boxno (local index) using the CSR format. function signature:
// [in ] gid gid[n] is the id for variable n at (i,j,k)
// [out] ncols # of non-zero columns in this row.
// [out] cols array of indices of columns with a non-zero matrix element in this row.
// [out] mat array of (non-zero) matrix elements in this row.
(int boxno, int i, int j, int k, int n, Array4<HYPRE_Int const> const* gid, HYPRE_Int& ncols, HYPRE_Int* cols, HYPRE_Real* mat)
[in]a_verboseHYPRE verbosity (default 0)
[in]a_options_namespacenamespace to use when parsing runtime input parameters (default "hypre")

Member Function Documentation

◆ fill_global_id()

template<int MSS>
void amrex::HypreSolver< MSS >::fill_global_id

◆ fill_local_id()

template<int MSS>
template<class Marker >
std::enable_if_t< IsCallable< Marker, int, int, int, int, int >::value > amrex::HypreSolver< MSS >::fill_local_id ( Marker const &  marker)

◆ fill_matrix()

template<int MSS>
template<class Filler , std::enable_if_t< IsCallable< Filler, int, int, int, int, int, Array4< HYPRE_Int const > const *, HYPRE_Int &, HYPRE_Int *, HYPRE_Real * >::value, int > FOO>
void amrex::HypreSolver< MSS >::fill_matrix ( Filler const &  filler)

◆ get_solution()

template<int MSS>
template<class MF , std::enable_if_t< IsFabArray< MF >::value &&std::is_same_v< typename MF::value_type, HYPRE_Real >, int > FOO>
void amrex::HypreSolver< MSS >::get_solution ( Vector< MF * > const &  a_soln)

◆ getA()

template<int MSS>
HYPRE_IJMatrix amrex::HypreSolver< MSS >::getA ( ) const
inline

◆ getb()

template<int MSS>
HYPRE_IJVector amrex::HypreSolver< MSS >::getb ( ) const
inline

◆ getFinalResidualNorm()

template<int MSS>
HYPRE_Real amrex::HypreSolver< MSS >::getFinalResidualNorm ( ) const
inline

◆ getNumIters()

template<int MSS>
int amrex::HypreSolver< MSS >::getNumIters ( ) const
inline

◆ getx()

template<int MSS>
HYPRE_IJVector amrex::HypreSolver< MSS >::getx ( ) const
inline

◆ load_vectors()

template<int MSS>
template<class MF , std::enable_if_t< IsFabArray< MF >::value &&std::is_same_v< typename MF::value_type, HYPRE_Real >, int > FOO>
void amrex::HypreSolver< MSS >::load_vectors ( Vector< MF * > const &  a_soln,
Vector< MF const * > const &  a_rhs 
)

◆ solve()

template<int MSS>
template<class MF , std::enable_if_t< IsFabArray< MF >::value &&std::is_same_v< typename MF::value_type, HYPRE_Real >, int > FOO>
void amrex::HypreSolver< MSS >::solve ( Vector< MF * > const &  a_soln,
Vector< MF const * > const &  a_rhs,
HYPRE_Real  rel_tol,
HYPRE_Real  abs_tol,
int  max_iter 
)

Solve Ax=b

Parameters
[out]a_solnthe solution (i.e. x)
[in]a_rhsthe right-hand side (i.e. b)
[in]rel_tolrelative convergence tolerance
[in]abs_tolabsolute convergence tolerance
[in]max_itermaximum number of iterations

Member Data Documentation

◆ m_A

template<int MSS>
HYPRE_IJMatrix amrex::HypreSolver< MSS >::m_A = nullptr
private

◆ m_b

template<int MSS>
HYPRE_IJVector amrex::HypreSolver< MSS >::m_b = nullptr
private

◆ m_cell_offset

template<int MSS>
LayoutData<Gpu::DeviceVector<int> > amrex::HypreSolver< MSS >::m_cell_offset
private

◆ m_comm

template<int MSS>
MPI_Comm amrex::HypreSolver< MSS >::m_comm = MPI_COMM_NULL
private

◆ m_dmap

template<int MSS>
DistributionMapping amrex::HypreSolver< MSS >::m_dmap
private

◆ m_geom

template<int MSS>
Geometry amrex::HypreSolver< MSS >::m_geom
private

◆ m_global_id

template<int MSS>
Vector<FabArray<BaseFab<HYPRE_Int> > > amrex::HypreSolver< MSS >::m_global_id
private

◆ m_global_id_vec

template<int MSS>
LayoutData<Gpu::DeviceVector<HYPRE_Int> > amrex::HypreSolver< MSS >::m_global_id_vec
private

◆ m_grids

template<int MSS>
Vector<BoxArray> amrex::HypreSolver< MSS >::m_grids
private

◆ m_hypre_ij

template<int MSS>
std::unique_ptr<HypreIJIface> amrex::HypreSolver< MSS >::m_hypre_ij
private

◆ m_id_offset

template<int MSS>
Vector<LayoutData<HYPRE_Int> > amrex::HypreSolver< MSS >::m_id_offset
private

◆ m_index_type

template<int MSS>
Vector<IndexType> amrex::HypreSolver< MSS >::m_index_type
private

◆ m_local_id

template<int MSS>
Vector<iMultiFab> amrex::HypreSolver< MSS >::m_local_id
private

◆ m_nghost

template<int MSS>
IntVect amrex::HypreSolver< MSS >::m_nghost
private

◆ m_nrows

template<int MSS>
LayoutData<HYPRE_Int> amrex::HypreSolver< MSS >::m_nrows
private

◆ m_nrows_grid

template<int MSS>
Vector<LayoutData<HYPRE_Int> > amrex::HypreSolver< MSS >::m_nrows_grid
private

◆ m_nrows_proc

template<int MSS>
HYPRE_Int amrex::HypreSolver< MSS >::m_nrows_proc
private

◆ m_nvars

template<int MSS>
int amrex::HypreSolver< MSS >::m_nvars
private

◆ m_options_namespace

template<int MSS>
std::string amrex::HypreSolver< MSS >::m_options_namespace
private

◆ m_owner_mask

template<int MSS>
Vector<std::unique_ptr<iMultiFab> > amrex::HypreSolver< MSS >::m_owner_mask
private

◆ m_verbose

template<int MSS>
int amrex::HypreSolver< MSS >::m_verbose
private

◆ m_x

template<int MSS>
HYPRE_IJVector amrex::HypreSolver< MSS >::m_x = nullptr
private

The documentation for this class was generated from the following file: