Block-Structured AMR Software Framework
amrex::GMRESMLMGT< MF > Class Template Reference

Solve using GMRES with multigrid as preconditioner. More...

#include <AMReX_GMRES_MLMG.H>

Public Types

using MG = MLMGT< MF >
 
using RT = typename MG::RT
 
using GM = GMRES< MF, GMRESMLMGT< MF > >
 

Public Member Functions

 GMRESMLMGT (MG &mlmg)
 
void solve (MF &a_sol, MF const &a_rhs, RT a_tol_rel, RT a_tol_abs)
 Solve the linear system. More...
 
void setVerbose (int v)
 Sets verbosity. More...
 
void setMaxIters (int niters)
 Sets the max number of iterations. More...
 
int getNumIters () const
 Gets the number of iterations. More...
 
RT getResidualNorm () const
 Gets the 2-norm of the residual. More...
 
GMgetGMRES ()
 Get the GMRES object. More...
 
void setPropertyOfZero (bool b)
 Set MLMG's multiplicative property of zero. More...
 
MF makeVecRHS () const
 Make MultiFab without ghost cells. More...
 
MF makeVecLHS () const
 Make MultiFab with ghost cells and set ghost cells to zero. More...
 
RT norm2 (MF const &mf) const
 
RT dotProduct (MF const &mf1, MF const &mf2) const
 
void apply (MF &lhs, MF const &rhs) const
 lhs = L(rhs) More...
 
void precond (MF &lhs, MF const &rhs) const
 
bool usePrecond (bool new_flag)
 Control whether or not to use MLMG as preconditioner. More...
 
void setPrecondNumIters (int precond_niters)
 Set the number of MLMG preconditioner iterations per GMRES iteration. More...
 

Static Public Member Functions

static void scale (MF &mf, RT scale_factor)
 
static void setToZero (MF &lhs)
 lhs = 0 More...
 
static void assign (MF &lhs, MF const &rhs)
 lhs = rhs More...
 
static void increment (MF &lhs, MF const &rhs, RT a)
 lhs += a*rhs More...
 
static void linComb (MF &lhs, RT a, MF const &rhs_a, RT b, MF const &rhs_b)
 lhs = a*rhs_a + b*rhs_b More...
 

Private Attributes

GM m_gmres
 
MGm_mlmg
 
MLLinOpT< MF > * m_linop
 
bool m_use_precond = true
 
bool m_prop_zero = false
 
int m_precond_niters = 1
 

Detailed Description

template<typename MF>
class amrex::GMRESMLMGT< MF >

Solve using GMRES with multigrid as preconditioner.

The linear system to solve is provided by MLMG, which is also being used as the preconditioner.

Member Typedef Documentation

◆ GM

template<typename MF >
using amrex::GMRESMLMGT< MF >::GM = GMRES<MF,GMRESMLMGT<MF> >

◆ MG

template<typename MF >
using amrex::GMRESMLMGT< MF >::MG = MLMGT<MF>

◆ RT

template<typename MF >
using amrex::GMRESMLMGT< MF >::RT = typename MG::RT

Constructor & Destructor Documentation

◆ GMRESMLMGT()

template<typename MF >
amrex::GMRESMLMGT< MF >::GMRESMLMGT ( MG mlmg)
explicit

Member Function Documentation

◆ apply()

template<typename MF >
void amrex::GMRESMLMGT< MF >::apply ( MF &  lhs,
MF const &  rhs 
) const

lhs = L(rhs)

◆ assign()

template<typename MF >
void amrex::GMRESMLMGT< MF >::assign ( MF &  lhs,
MF const &  rhs 
)
static

lhs = rhs

◆ dotProduct()

template<typename MF >
auto amrex::GMRESMLMGT< MF >::dotProduct ( MF const &  mf1,
MF const &  mf2 
) const

◆ getGMRES()

template<typename MF >
GM& amrex::GMRESMLMGT< MF >::getGMRES ( )
inline

Get the GMRES object.

◆ getNumIters()

template<typename MF >
int amrex::GMRESMLMGT< MF >::getNumIters ( ) const
inline

Gets the number of iterations.

◆ getResidualNorm()

template<typename MF >
RT amrex::GMRESMLMGT< MF >::getResidualNorm ( ) const
inline

Gets the 2-norm of the residual.

◆ increment()

template<typename MF >
void amrex::GMRESMLMGT< MF >::increment ( MF &  lhs,
MF const &  rhs,
RT  a 
)
static

lhs += a*rhs

◆ linComb()

template<typename MF >
void amrex::GMRESMLMGT< MF >::linComb ( MF &  lhs,
RT  a,
MF const &  rhs_a,
RT  b,
MF const &  rhs_b 
)
static

lhs = a*rhs_a + b*rhs_b

◆ makeVecLHS()

template<typename MF >
auto amrex::GMRESMLMGT< MF >::makeVecLHS

Make MultiFab with ghost cells and set ghost cells to zero.

◆ makeVecRHS()

template<typename MF >
auto amrex::GMRESMLMGT< MF >::makeVecRHS

Make MultiFab without ghost cells.

◆ norm2()

template<typename MF >
auto amrex::GMRESMLMGT< MF >::norm2 ( MF const &  mf) const

◆ precond()

template<typename MF >
void amrex::GMRESMLMGT< MF >::precond ( MF &  lhs,
MF const &  rhs 
) const

◆ scale()

template<typename MF >
void amrex::GMRESMLMGT< MF >::scale ( MF &  mf,
RT  scale_factor 
)
static

◆ setMaxIters()

template<typename MF >
void amrex::GMRESMLMGT< MF >::setMaxIters ( int  niters)
inline

Sets the max number of iterations.

◆ setPrecondNumIters()

template<typename MF >
void amrex::GMRESMLMGT< MF >::setPrecondNumIters ( int  precond_niters)
inline

Set the number of MLMG preconditioner iterations per GMRES iteration.

◆ setPropertyOfZero()

template<typename MF >
void amrex::GMRESMLMGT< MF >::setPropertyOfZero ( bool  b)
inline

Set MLMG's multiplicative property of zero.

This should NOT be called unless MLMG has the multiplicative property of zero. MLMG is not a matrix, and usually does not have the properties of a matrix such as the multiplicative property of zero (i.e., M*0=0) because how domain boundary conditions are handled. However, if MLMG has the property of zero, calling this function with true can have a small performance benefit.

◆ setToZero()

template<typename MF >
void amrex::GMRESMLMGT< MF >::setToZero ( MF &  lhs)
static

lhs = 0

◆ setVerbose()

template<typename MF >
void amrex::GMRESMLMGT< MF >::setVerbose ( int  v)
inline

Sets verbosity.

◆ solve()

template<typename MF >
void amrex::GMRESMLMGT< MF >::solve ( MF &  a_sol,
MF const &  a_rhs,
RT  a_tol_rel,
RT  a_tol_abs 
)

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.

◆ usePrecond()

template<typename MF >
bool amrex::GMRESMLMGT< MF >::usePrecond ( bool  new_flag)
inline

Control whether or not to use MLMG as preconditioner.

Member Data Documentation

◆ m_gmres

template<typename MF >
GM amrex::GMRESMLMGT< MF >::m_gmres
private

◆ m_linop

template<typename MF >
MLLinOpT<MF>* amrex::GMRESMLMGT< MF >::m_linop
private

◆ m_mlmg

template<typename MF >
MG* amrex::GMRESMLMGT< MF >::m_mlmg
private

◆ m_precond_niters

template<typename MF >
int amrex::GMRESMLMGT< MF >::m_precond_niters = 1
private

◆ m_prop_zero

template<typename MF >
bool amrex::GMRESMLMGT< MF >::m_prop_zero = false
private

◆ m_use_precond

template<typename MF >
bool amrex::GMRESMLMGT< MF >::m_use_precond = true
private

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