Block-Structured AMR Software Framework
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)
 Wrap an existing MLMG hierarchy as a GMRES operator/preconditioner.
 
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)
 Solve using pointer-to-MultiFab arrays (utility for existing MLMG data).
 
void setVerbose (int v)
 Set verbosity level v for both GMRES and MLMG components.
 
void setMaxIters (int niters)
 Cap the number of GMRES iterations executed per solve to niters.
 
int getNumIters () const
 Number of iterations executed by the last solve().
 
RT getResidualNorm () const
 Final residual 2-norm produced by the last solve().
 
GMgetGMRES ()
 Direct access to the underlying GMRES driver for fine-grained tuning.
 
VEC makeVecRHS () const
 Return a MultiFab vector without ghost cells for RHS storage.
 
VEC makeVecLHS () const
 
RT norm2 (VEC const &mf) const
 Return the 2-norm of mf using the multigrid operator's weighting.
 
RT dotProduct (VEC const &mf1, VEC const &mf2) const
 Return the preconditioner-aware dot product between mf1 and mf2.
 
void apply (VEC &lhs, VEC const &rhs) const
 Apply the operator supplied to the constructor.
 
void precond (VEC &lhs, VEC const &rhs) const
 Apply the MLMG-based preconditioner: $\text{lhs} = P^{-1}(\text{rhs})$.
 
bool usePrecond (bool new_flag)
 Enable or disable MLMG as a preconditioner.
 
void setPrecondNumIters (int precond_niters)
 Adjust how many MLMG smoothing iterations run inside each GMRES iteration.
 

Static Public Member Functions

static void scale (VEC &mf, RT scale_factor)
 Scale each component of mf by scale_factor.
 
static void setToZero (VEC &lhs)
 Reset lhs to zero on every AMR level.
 
static void assign (VEC &lhs, VEC const &rhs)
 Copy rhs into lhs level by level.
 
static void increment (VEC &lhs, VEC const &rhs, RT a)
 Add a scaled vector: $ \text{lhs} \mathrel{+}= a\,\text{rhs}$.
 
static void linComb (VEC &lhs, RT a, VEC const &rhs_a, RT b, VEC const &rhs_b)
 Form a linear combination: $ \text{lhs} = a\,\text{rhs}_a + b\,\text{rhs}_b$.
 

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

Wrap an existing MLMG hierarchy as a GMRES operator/preconditioner.

Parameters
mlmgMultigrid object that provides apply/precond hooks.

Member Function Documentation

◆ apply()

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

Apply the operator supplied to the constructor.

Parameters
lhsDestination vector receiving $L(\text{rhs})$.
rhsSource vector.

◆ assign()

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

Copy rhs into lhs level by level.

Parameters
lhsDestination vector updated in place.
rhsSource vector.

◆ dotProduct()

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

Return the preconditioner-aware dot product between mf1 and mf2.

◆ getGMRES()

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

Direct access to the underlying GMRES driver for fine-grained tuning.

◆ getNumIters()

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

Number of iterations executed by the last solve().

◆ getResidualNorm()

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

Final residual 2-norm produced by the last solve().

◆ increment()

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

Add a scaled vector: $ \text{lhs} \mathrel{+}= a\,\text{rhs}$.

Parameters
lhsAccumulation target.
rhsSource vector.
aScalar multiplier applied to 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

Form a linear combination: $ \text{lhs} = a\,\text{rhs}_a + b\,\text{rhs}_b$.

Parameters
lhsDestination vector written in place.
aScalar applied to rhs_a.
rhs_aFirst source vector.
bScalar applied to rhs_b.
rhs_bSecond source vector.

◆ makeVecLHS()

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

Return a MultiFab vector with one layer of ghost cells for LHS storage. Also set the data in the ghost cells to zero.

◆ makeVecRHS()

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

Return a MultiFab vector without ghost cells for RHS storage.

◆ norm2()

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

Return the 2-norm of mf using the multigrid operator's weighting.

◆ precond()

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

Apply the MLMG-based preconditioner: $\text{lhs} = P^{-1}(\text{rhs})$.

Parameters
lhsDestination vector, zeroed before solving when the preconditioner is enabled.
rhsSource vector passed to the preconditioner.

◆ scale()

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

Scale each component of mf by scale_factor.

◆ setMaxIters()

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

Cap the number of GMRES iterations executed per solve to niters.

◆ setPrecondNumIters()

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

Adjust how many MLMG smoothing iterations run inside each GMRES iteration.

Parameters
precond_nitersNumber of MLMG iterations per GMRES iteration.

◆ setToZero()

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

Reset lhs to zero on every AMR level.

Parameters
lhsPer-level vector that will be zeroed in place.

◆ setVerbose()

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

Set verbosity level v for both GMRES and MLMG components.

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

Solve using pointer-to-MultiFab arrays (utility for existing MLMG data).

Parameters
a_solPer-level solution MultiFabs (single component).
a_rhsPer-level RHS MultiFabs (single component).
a_tol_relRelative tolerance.
a_tol_absAbsolute tolerance.

◆ usePrecond()

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

Enable or disable MLMG as a preconditioner.

Parameters
new_flagSet true to enable, false to disable.
Returns
Previous flag value.

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