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< IndexType > | m_index_type |
IntVect | m_nghost |
Geometry | m_geom |
Vector< BoxArray > | m_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< iMultiFab > | m_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< HypreIJIface > | m_hypre_ij |
HYPRE_IJMatrix | m_A = nullptr |
HYPRE_IJVector | m_b = nullptr |
HYPRE_IJVector | m_x = nullptr |
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
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
[in] | a_index_type | nodality of the data |
[in] | a_nghost | number of ghosts cells |
[in] | a_geom | Geometry of problem domain |
[in] | a_grids | BoxArray of computational grids |
[in] | a_dmap | DistributionMapping |
[in] | a_marker | functor 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_filler | functor 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_verbose | HYPRE verbosity (default 0) |
[in] | a_options_namespace | namespace to use when parsing runtime input parameters (default "hypre") |
void amrex::HypreSolver< MSS >::fill_global_id |
std::enable_if_t< IsCallable< Marker, int, int, int, int, int >::value > amrex::HypreSolver< MSS >::fill_local_id | ( | Marker const & | marker | ) |
void amrex::HypreSolver< MSS >::fill_matrix | ( | Filler const & | filler | ) |
void amrex::HypreSolver< MSS >::get_solution | ( | Vector< MF * > const & | a_soln | ) |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
void amrex::HypreSolver< MSS >::load_vectors | ( | Vector< MF * > const & | a_soln, |
Vector< MF const * > const & | a_rhs | ||
) |
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
[out] | a_soln | the solution (i.e. x) |
[in] | a_rhs | the right-hand side (i.e. b) |
[in] | rel_tol | relative convergence tolerance |
[in] | abs_tol | absolute convergence tolerance |
[in] | max_iter | maximum number of iterations |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |