![]() |
Block-Structured AMR Software Framework
|
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(). | |
| GM & | getGMRES () |
| 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: ![]() | |
| 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: ![]() | |
| static void | linComb (VEC &lhs, RT a, VEC const &rhs_a, RT b, VEC const &rhs_b) |
Form a linear combination: ![]() | |
Solve using GMRES with multigrid as preconditioner.
The linear system to solve is provided by MLMG, which is also being used as the preconditioner.
| using amrex::GMRESMLMGT< MF >::GM = GMRES<Vector<MF>,GMRESMLMGT<MF> > |
| using amrex::GMRESMLMGT< MF >::MG = MLMGT<MF> |
| using amrex::GMRESMLMGT< MF >::RT = typename MG::RT |
| using amrex::GMRESMLMGT< MF >::VEC = Vector<MF> |
|
explicit |
Wrap an existing MLMG hierarchy as a GMRES operator/preconditioner.
| mlmg | Multigrid object that provides apply/precond hooks. |
| void amrex::GMRESMLMGT< MF >::apply | ( | VEC & | lhs, |
| VEC const & | rhs | ||
| ) | const |
Apply the operator supplied to the constructor.
| lhs | Destination vector receiving ![]() |
| rhs | Source vector. |
|
static |
Copy rhs into lhs level by level.
| lhs | Destination vector updated in place. |
| rhs | Source vector. |
| auto amrex::GMRESMLMGT< MF >::dotProduct | ( | VEC const & | mf1, |
| VEC const & | mf2 | ||
| ) | const |
Return the preconditioner-aware dot product between mf1 and mf2.
|
inline |
Direct access to the underlying GMRES driver for fine-grained tuning.
|
inline |
Number of iterations executed by the last solve().
|
inline |
Final residual 2-norm produced by the last solve().
|
static |
Add a scaled vector: 
| lhs | Accumulation target. |
| rhs | Source vector. |
| a | Scalar multiplier applied to rhs. |
|
static |
Form a linear combination: 
| lhs | Destination vector written in place. |
| a | Scalar applied to rhs_a. |
| rhs_a | First source vector. |
| b | Scalar applied to rhs_b. |
| rhs_b | Second source vector. |
| 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.
| auto amrex::GMRESMLMGT< MF >::makeVecRHS | ( | ) | const |
Return a MultiFab vector without ghost cells for RHS storage.
| auto amrex::GMRESMLMGT< MF >::norm2 | ( | VEC const & | mf | ) | const |
Return the 2-norm of mf using the multigrid operator's weighting.
| void amrex::GMRESMLMGT< MF >::precond | ( | VEC & | lhs, |
| VEC const & | rhs | ||
| ) | const |
Apply the MLMG-based preconditioner: 
| lhs | Destination vector, zeroed before solving when the preconditioner is enabled. |
| rhs | Source vector passed to the preconditioner. |
|
static |
Scale each component of mf by scale_factor.
|
inline |
Cap the number of GMRES iterations executed per solve to niters.
|
inline |
|
static |
Reset lhs to zero on every AMR level.
| lhs | Per-level vector that will be zeroed in place. |
|
inline |
Set verbosity level v for both GMRES and MLMG components.
| void amrex::GMRESMLMGT< MF >::solve | ( | MF & | a_sol, |
| MF const & | a_rhs, | ||
| RT | a_tol_rel, | ||
| RT | a_tol_abs | ||
| ) |
Solve the linear system.
| a_sol | unknowns, i.e., x in A x = b. |
| a_rhs | RHS, i.e., b in A x = b. |
| a_tol_rel | relative tolerance. |
| a_tol_abs | absolute tolerance. |
| 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).
| a_sol | Per-level solution MultiFabs (single component). |
| a_rhs | Per-level RHS MultiFabs (single component). |
| a_tol_rel | Relative tolerance. |
| a_tol_abs | Absolute tolerance. |
|
inline |
Enable or disable MLMG as a preconditioner.
| new_flag | Set true to enable, false to disable. |