Block-Structured AMR Software Framework
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
Loading...
Searching...
No Matches
amrex::GMRESMLMGT< MF > Class Template Reference

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

#include <AMReX_GMRES_MLMG.H>

Public Types

using VEC = Vector< MF >
 
using MG = MLMGT< MF >
 
using RT = typename MG::RT
 
using GM = GMRES< Vector< 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.
 
void solve (Vector< MF * > const &a_sol, Vector< MF const * > const &a_rhs, RT a_tol_rel, RT a_tol_abs)
 
void setVerbose (int v)
 Sets verbosity.
 
void setMaxIters (int niters)
 Sets the max number of iterations.
 
int getNumIters () const
 Gets the number of iterations.
 
RT getResidualNorm () const
 Gets the 2-norm of the residual.
 
GMgetGMRES ()
 Get the GMRES object.
 
VEC makeVecRHS () const
 Make MultiFab without ghost cells.
 
VEC makeVecLHS () const
 Make MultiFab with ghost cells and set ghost cells to zero.
 
RT norm2 (VEC const &mf) const
 
RT dotProduct (VEC const &mf1, VEC const &mf2) const
 
void apply (VEC &lhs, VEC const &rhs) const
 lhs = L(rhs)
 
void precond (VEC &lhs, VEC const &rhs) const
 
bool usePrecond (bool new_flag)
 Control whether or not to use MLMG as preconditioner.
 
void setPrecondNumIters (int precond_niters)
 Set the number of MLMG preconditioner iterations per GMRES iteration.
 

Static Public Member Functions

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

Private Attributes

GM m_gmres
 
MGm_mlmg
 
MLLinOpT< MF > * m_linop
 
int m_verbose = 0
 
int m_nlevels = 0
 
bool m_use_precond = true
 
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<Vector<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

◆ VEC

template<typename MF >
using amrex::GMRESMLMGT< MF >::VEC = Vector<MF>

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 ( VEC lhs,
VEC const &  rhs 
) const

lhs = L(rhs)

◆ assign()

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

lhs = rhs

◆ dotProduct()

template<typename MF >
auto amrex::GMRESMLMGT< MF >::dotProduct ( VEC const &  mf1,
VEC 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 ( VEC lhs,
VEC const &  rhs,
RT  a 
)
static

lhs += a*rhs

◆ linComb()

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

lhs = a*rhs_a + b*rhs_b

◆ makeVecLHS()

template<typename MF >
auto amrex::GMRESMLMGT< MF >::makeVecLHS ( ) const

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

◆ makeVecRHS()

template<typename MF >
auto amrex::GMRESMLMGT< MF >::makeVecRHS ( ) const

Make MultiFab without ghost cells.

◆ norm2()

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

◆ precond()

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

◆ scale()

template<typename MF >
void amrex::GMRESMLMGT< MF >::scale ( VEC 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.

◆ setToZero()

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

lhs = 0

◆ setVerbose()

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

Sets verbosity.

◆ solve() [1/2]

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.

◆ solve() [2/2]

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

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

template<typename MF >
int amrex::GMRESMLMGT< MF >::m_nlevels = 0
private

◆ m_precond_niters

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

◆ m_use_precond

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

◆ m_verbose

template<typename MF >
int amrex::GMRESMLMGT< MF >::m_verbose = 0
private

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