Block-Structured AMR Software Framework
amrex::GMRES< V, M > Class Template Reference

GMRES. More...

#include <AMReX_GMRES.H>

Public Types

using RT = typename M::RT
 

Public Member Functions

 GMRES ()
 
void define (M &linop)
 
void solve (V &a_sol, V const &a_rhs, RT a_tol_rel, RT a_tol_abs, int a_its=-1)
 Solve the linear system. More...
 
void setVerbose (int v)
 Sets verbosity. More...
 
void setRestartLength (int rl)
 Sets restart length. The default is 30. More...
 
void setMaxIters (int niters)
 Sets the max number of iterations. More...
 
int getNumIters () const
 Gets the number of iterations. More...
 
int getStatus () const
 Gets the solver status. More...
 
RT getResidualNorm () const
 Gets the 2-norm of the residual. More...
 

Private Member Functions

void clear ()
 
void allocate_scratch ()
 
void cycle (V &a_xx, int &a_status, int &a_itcount, RT &a_rnorm0)
 
void build_solution (V &a_xx, int it)
 
void compute_residual (V &a_rr, V const &a_xx, V const &a_bb)
 
bool converged (RT r0, RT r) const
 
void gram_schmidt_orthogonalization (int it)
 
void update_hessenberg (int it, bool happyend, RT &res)
 

Private Attributes

int m_verbose = 0
 
int m_maxiter = 2000
 
int m_its = 0
 
int m_status = -1
 
int m_restrtlen = 30
 
RT m_res = std::numeric_limits<RT>::max()
 
RT m_rtol = RT(0)
 
RT m_atol = RT(0)
 
Vector< RTm_hh_1d
 
Vector< RTm_hes_1d
 
Table2D< RTm_hh
 
Table2D< RTm_hes
 
Vector< RTm_grs
 
Vector< RTm_cc
 
Vector< RTm_ss
 
std::unique_ptr< V > m_v_tmp_rhs
 
std::unique_ptr< V > m_v_tmp_lhs
 
Vector< V > m_vv
 
M * m_linop = nullptr
 

Detailed Description

template<typename V, typename M>
class amrex::GMRES< V, M >

GMRES.

This class implements the GMRES algorithm. The template parameter V is for a linear algebra vector class. For example, it could be amrex::MultiFab. The other template parameter M is for a linear operator class with a number of required member functions. Note that conceptually M contains a matrix. However, it does not mean it needs to have a member variable storing the matrix, because GMRES only needs the matrix vector product, not the matrix itself.

Template Parameters
Vlinear algebra vector. It must be default constructible, move constructible, and move assignable.
Mlinear operator. A list of required member functions for M is shown below. Here RT (typename M::RT) is either double or float. Examples of implementation for V being a single-component MultiFab are also given below.
  • void apply(V& Ax, V const& x)
    Ax = A(x), where A is the linear operator performing matrix vector product. Here x is made with the makeVecLHS function. Therefore, it may have ghost cells and it's safe to cast it with const_cast<V&>(x) and do ghost cell exchange, if necessary.
  • void assign(V& lhs, V const& rhs)
    lhs = rhs. For example, MultiFab::Copy(lhs,rhs,0,0,1,0).
  • RT dotProduct(V const& v1, V const& v2)
    returns v1 * v2. For example, MultiFab::Dot(v1,0,v2,0,1,0).
  • void increment(V& lhs, V const& rhs, RT a)
    lhs += a * rhs. For example, MultiFab::Saxpy(lhs,a,rhs,0,0,1,0).
  • void linComb(V& lhs, RT a, V const& rhs_a, RT b, V const& rhs_b)
    lhs = a * rhs_a + b * rhs_b. For example, MultiFab::LinComb(lhs,a,rhs_a,0,b,rhs_b,0,0,1,0).
  • V makeVecRHS()
    returns a V object that is suitable as RHS in M x = b. The reason we distinguish between LHS and RHS is M might need the distinction for efficiency. For example, if V is MultiFab, we might need the x in the LHS of M x = b to have ghost cells for efficiency, whereas no ghost cells are needed for the RHS (i.e., b). An example of the implementation might be return MultiFab(grids,dmap,1,0).
  • V makeVecLHS()
    returns a V object that is suitable as LHS in M x = b. See the description for makeVecRHS for more details. An example of the implementation might be return MultiFab(grids,dmap,1,1).
  • RT norm2(V const& v)
    returns the 2-norm of v. For example, return v.norm2().
  • void precond(V& lhs, V const& rhs)
    applies right-preconditioning, i.e., solve P(lhs) = rhs, where P is an approximation to A If there is no preconditioner, P = I and thus this function should do lhs = rhs.
  • void scale(V& v, RT fac)
    scales v by fac. For example, v.mult(fac).
  • void setToZero(V& v)
    v = 0. For example, v.setVal(0).

Member Typedef Documentation

◆ RT

template<typename V , typename M >
using amrex::GMRES< V, M >::RT = typename M::RT

Constructor & Destructor Documentation

◆ GMRES()

template<typename V , typename M >
amrex::GMRES< V, M >::GMRES

Member Function Documentation

◆ allocate_scratch()

template<typename V , typename M >
void amrex::GMRES< V, M >::allocate_scratch
private

◆ build_solution()

template<typename V , typename M >
void amrex::GMRES< V, M >::build_solution ( V &  a_xx,
int  it 
)
private

◆ clear()

template<typename V , typename M >
void amrex::GMRES< V, M >::clear
private

◆ compute_residual()

template<typename V , typename M >
void amrex::GMRES< V, M >::compute_residual ( V &  a_rr,
V const &  a_xx,
V const &  a_bb 
)
private

◆ converged()

template<typename V , typename M >
bool amrex::GMRES< V, M >::converged ( RT  r0,
RT  r 
) const
private

◆ cycle()

template<typename V , typename M >
void amrex::GMRES< V, M >::cycle ( V &  a_xx,
int a_status,
int a_itcount,
RT a_rnorm0 
)
private

◆ define()

template<typename V , typename M >
void amrex::GMRES< V, M >::define ( M &  linop)

Defines with a reference to M. It's the user's responsibility to keep the M object alive for GMRES to be functional. This function must be called before solve() can be called.

◆ getNumIters()

template<typename V , typename M >
int amrex::GMRES< V, M >::getNumIters ( ) const
inline

Gets the number of iterations.

◆ getResidualNorm()

template<typename V , typename M >
RT amrex::GMRES< V, M >::getResidualNorm ( ) const
inline

Gets the 2-norm of the residual.

◆ getStatus()

template<typename V , typename M >
int amrex::GMRES< V, M >::getStatus ( ) const
inline

Gets the solver status.

◆ gram_schmidt_orthogonalization()

template<typename V , typename M >
void amrex::GMRES< V, M >::gram_schmidt_orthogonalization ( int  it)
private

◆ setMaxIters()

template<typename V , typename M >
void amrex::GMRES< V, M >::setMaxIters ( int  niters)
inline

Sets the max number of iterations.

◆ setRestartLength()

template<typename V , typename M >
void amrex::GMRES< V, M >::setRestartLength ( int  rl)

Sets restart length. The default is 30.

◆ setVerbose()

template<typename V , typename M >
void amrex::GMRES< V, M >::setVerbose ( int  v)
inline

Sets verbosity.

◆ solve()

template<typename V , typename M >
void amrex::GMRES< V, M >::solve ( V &  a_sol,
V const &  a_rhs,
RT  a_tol_rel,
RT  a_tol_abs,
int  a_its = -1 
)

Solve the linear system.

Parameters
a_solunknowns, i.e., x in A x = b.
a_rhsRHS, i.e., b in A x = b.
a_tol_relrelative tolerance.
a_tol_absabsolute tolerance. \patam a_its optional argument specifying the maximum number of iterations.

◆ update_hessenberg()

template<typename V , typename M >
void amrex::GMRES< V, M >::update_hessenberg ( int  it,
bool  happyend,
RT res 
)
private

Member Data Documentation

◆ m_atol

template<typename V , typename M >
RT amrex::GMRES< V, M >::m_atol = RT(0)
private

◆ m_cc

template<typename V , typename M >
Vector<RT> amrex::GMRES< V, M >::m_cc
private

◆ m_grs

template<typename V , typename M >
Vector<RT> amrex::GMRES< V, M >::m_grs
private

◆ m_hes

template<typename V , typename M >
Table2D<RT> amrex::GMRES< V, M >::m_hes
private

◆ m_hes_1d

template<typename V , typename M >
Vector<RT> amrex::GMRES< V, M >::m_hes_1d
private

◆ m_hh

template<typename V , typename M >
Table2D<RT> amrex::GMRES< V, M >::m_hh
private

◆ m_hh_1d

template<typename V , typename M >
Vector<RT> amrex::GMRES< V, M >::m_hh_1d
private

◆ m_its

template<typename V , typename M >
int amrex::GMRES< V, M >::m_its = 0
private

◆ m_linop

template<typename V , typename M >
M* amrex::GMRES< V, M >::m_linop = nullptr
private

◆ m_maxiter

template<typename V , typename M >
int amrex::GMRES< V, M >::m_maxiter = 2000
private

◆ m_res

template<typename V , typename M >
RT amrex::GMRES< V, M >::m_res = std::numeric_limits<RT>::max()
private

◆ m_restrtlen

template<typename V , typename M >
int amrex::GMRES< V, M >::m_restrtlen = 30
private

◆ m_rtol

template<typename V , typename M >
RT amrex::GMRES< V, M >::m_rtol = RT(0)
private

◆ m_ss

template<typename V , typename M >
Vector<RT> amrex::GMRES< V, M >::m_ss
private

◆ m_status

template<typename V , typename M >
int amrex::GMRES< V, M >::m_status = -1
private

◆ m_v_tmp_lhs

template<typename V , typename M >
std::unique_ptr<V> amrex::GMRES< V, M >::m_v_tmp_lhs
private

◆ m_v_tmp_rhs

template<typename V , typename M >
std::unique_ptr<V> amrex::GMRES< V, M >::m_v_tmp_rhs
private

◆ m_verbose

template<typename V , typename M >
int amrex::GMRES< V, M >::m_verbose = 0
private

◆ m_vv

template<typename V , typename M >
Vector<V> amrex::GMRES< V, M >::m_vv
private

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