#include <AMReX_MLMG.H>
Classes | |
| class | error |
Public Types | |
| enum class | CFStrategy : int { none , ghostnodes } |
| using | MFType = MF |
| using | FAB = typename MLLinOpT< MF >::FAB |
| using | RT = typename MLLinOpT< MF >::RT |
| using | BCMode = typename MLLinOpT< MF >::BCMode |
| using | Location = typename MLLinOpT< MF >::Location |
| using | BottomSolver = amrex::BottomSolver |
Public Member Functions | |
| MLMGT (MLLinOpT< MF > &a_lp) | |
| ~MLMGT () | |
| MLMGT (MLMGT< MF > const &)=delete | |
| MLMGT (MLMGT< MF > &&)=delete | |
| MLMGT< MF > & | operator= (MLMGT< MF > const &)=delete |
| MLMGT< MF > & | operator= (MLMGT< MF > &&)=delete |
| template<typename AMF > | |
| RT | solve (const Vector< AMF * > &a_sol, const Vector< AMF const * > &a_rhs, RT a_tol_rel, RT a_tol_abs, const char *checkpoint_file=nullptr) |
| template<typename AMF > | |
| RT | solve (std::initializer_list< AMF * > a_sol, std::initializer_list< AMF const * > a_rhs, RT a_tol_rel, RT a_tol_abs, const char *checkpoint_file=nullptr) |
| RT | precond (Vector< MF * > const &a_sol, Vector< MF const * > const &a_rhs, RT a_tol_rel, RT a_tol_abs) |
| template<typename AMF > | |
| void | getGradSolution (const Vector< Array< AMF *, 3 > > &a_grad_sol, Location a_loc=Location::FaceCenter) |
| template<typename AMF > | |
| void | getGradSolution (std::initializer_list< Array< AMF *, 3 > > a_grad_sol, Location a_loc=Location::FaceCenter) |
| template<typename AMF > | |
| void | getFluxes (const Vector< Array< AMF *, 3 > > &a_flux, Location a_loc=Location::FaceCenter) |
For (alpha * a - beta * (del dot b grad)) phi = rhs, flux means -b grad phi | |
| template<typename AMF > | |
| void | getFluxes (std::initializer_list< Array< AMF *, 3 > > a_flux, Location a_loc=Location::FaceCenter) |
| template<typename AMF > | |
| void | getFluxes (const Vector< Array< AMF *, 3 > > &a_flux, const Vector< AMF * > &a_sol, Location a_loc=Location::FaceCenter) |
| template<typename AMF > | |
| void | getFluxes (std::initializer_list< Array< AMF *, 3 > > a_flux, std::initializer_list< AMF * > a_sol, Location a_loc=Location::FaceCenter) |
| template<typename AMF > | |
| void | getFluxes (const Vector< AMF * > &a_flux, Location a_loc=Location::CellCenter) |
| template<typename AMF > | |
| void | getFluxes (std::initializer_list< AMF * > a_flux, Location a_loc=Location::CellCenter) |
| template<typename AMF > | |
| void | getFluxes (const Vector< AMF * > &a_flux, const Vector< AMF * > &a_sol, Location a_loc=Location::CellCenter) |
| template<typename AMF > | |
| void | getFluxes (std::initializer_list< AMF * > a_flux, std::initializer_list< AMF * > a_sol, Location a_loc=Location::CellCenter) |
| void | compResidual (const Vector< MF * > &a_res, const Vector< MF * > &a_sol, const Vector< MF const * > &a_rhs) |
| void | getEBFluxes (const Vector< MF * > &a_eb_flux) |
| void | getEBFluxes (const Vector< MF * > &a_eb_flux, const Vector< MF * > &a_sol) |
| void | apply (const Vector< MF * > &out, const Vector< MF * > &in) |
out = L(in). Note that, if no actual solve is needed, one could turn off multigrid coarsening by constructing a MLLinOp object with an appropriate LPInfo object (e.g., with LPInfo().setMaxCoarseningLevel(0)). | |
| void | applyPrecond (const Vector< MF * > &out, const Vector< MF * > &in) |
| out = L(in) as a preconditioner | |
| int | getVerbose () const |
| int | getBottomVerbose () const |
| void | incPrintIdentation () |
| void | decPrintIdentation () |
| void | setThrowException (bool t) noexcept |
| void | setVerbose (int v) noexcept |
| void | setMaxIter (int n) noexcept |
| void | setMaxFmgIter (int n) noexcept |
| void | setFixedIter (int nit) noexcept |
| void | setPrecondIter (int nit) noexcept |
| void | setPreSmooth (int n) noexcept |
| void | setPostSmooth (int n) noexcept |
| void | setFinalSmooth (int n) noexcept |
| void | setBottomSmooth (int n) noexcept |
| void | setBottomSolver (BottomSolver s) noexcept |
| BottomSolver | getBottomSolver () const noexcept |
| void | setCFStrategy (CFStrategy a_cf_strategy) noexcept |
| void | setBottomVerbose (int v) noexcept |
| void | setBottomMaxIter (int n) noexcept |
| void | setBottomTolerance (RT t) noexcept |
| void | setBottomToleranceAbs (RT t) noexcept |
| RT | getBottomToleranceAbs () const noexcept |
| void | setAlwaysUseBNorm (int flag) noexcept |
| void | setConvergenceNormType (MLMGNormType norm) noexcept |
| void | setFinalFillBC (int flag) noexcept |
| int | numAMRLevels () const noexcept |
| void | setNSolve (int flag) noexcept |
| void | setNSolveGridSize (int s) noexcept |
| void | setNoGpuSync (bool do_not_sync) noexcept |
| void | prepareForFluxes (Vector< MF const * > const &a_sol) |
| template<typename AMF > | |
| void | prepareForSolve (Vector< AMF * > const &a_sol, Vector< AMF const * > const &a_rhs) |
| void | prepareForNSolve () |
| void | prepareLinOp () |
| void | preparePrecond () |
| void | oneIter (int iter) |
| void | miniCycle (int amrlev) |
| void | mgVcycle (int amrlev, int mglev) |
| void | mgFcycle () |
| void | bottomSolve () |
| void | NSolve (MLMGT< MF > &a_solver, MF &a_sol, MF &a_rhs) |
| void | actualBottomSolve () |
| void | computeMLResidual (int amrlevmax) |
| void | computeResidual (int alev) |
| void | computeResWithCrseSolFineCor (int calev, int falev) |
| void | computeResWithCrseCorFineCor (int falev) |
| void | interpCorrection (int alev) |
| void | interpCorrection (int alev, int mglev) |
| void | addInterpCorrection (int alev, int mglev) |
| void | computeResOfCorrection (int amrlev, int mglev) |
| RT | ResNormInf (int alev, bool local=false) |
| RT | MLResNormInf (int alevmax, bool local=false) |
| RT | MLRhsNormInf (bool local=false) |
| void | makeSolvable () |
| void | makeSolvable (int amrlev, int mglev, MF &mf) |
| int | bottomSolveWithCG (MF &x, const MF &b, typename MLCGSolverT< MF >::Type type) |
| RT | getInitRHS () const noexcept |
| RT | getInitResidual () const noexcept |
| RT | getFinalResidual () const noexcept |
| Vector< RT > const & | getResidualHistory () const noexcept |
| int | getNumIters () const noexcept |
| Vector< int > const & | getNumCGIters () const noexcept |
| MLLinOpT< MF > & | getLinOp () |
| template<typename AMF > | |
| auto | solve (std::initializer_list< AMF * > a_sol, std::initializer_list< AMF const * > a_rhs, RT a_tol_rel, RT a_tol_abs, const char *checkpoint_file) -> RT |
| template<typename AMF > | |
| auto | solve (const Vector< AMF * > &a_sol, const Vector< AMF const * > &a_rhs, RT a_tol_rel, RT a_tol_abs, const char *checkpoint_file) -> RT |
Private Types | |
| enum | timer_types { solve_time =0 , iter_time , bottom_time , ntimers } |
Private Member Functions | |
| void | checkPoint (const Vector< MultiFab * > &a_sol, const Vector< MultiFab const * > &a_rhs, RT a_tol_rel, RT a_tol_abs, const char *a_file_name) const |
Private Attributes | |
| bool | precond_mode = false |
| bool | throw_exception = false |
| int | verbose = 1 |
| int | max_iters = 200 |
| int | do_fixed_number_of_iters = 0 |
| int | max_precond_iters = 1 |
| int | nu1 = 2 |
| pre | |
| int | nu2 = 2 |
| post | |
| int | nuf = 8 |
| when smoother is used as bottom solver | |
| int | nub = 0 |
| additional smoothing after bottom cg solver | |
| int | max_fmg_iters = 0 |
| BottomSolver | bottom_solver = BottomSolver::Default |
| CFStrategy | cf_strategy = CFStrategy::none |
| int | bottom_verbose = 0 |
| int | bottom_maxiter = 200 |
| RT | bottom_reltol = std::is_same<RT,double>() ? RT(1.e-4) : RT(1.e-3) |
| RT | bottom_abstol = RT(-1.0) |
| MLMGNormType | norm_type = MLMGNormType::greater |
| int | final_fill_bc = 0 |
| MLLinOpT< MF > & | linop |
| int | ncomp |
| int | namrlevs |
| int | finest_amr_lev |
| bool | linop_prepared = false |
| Long | solve_called = 0 |
| int | do_nsolve = false |
| N Solve. | |
| int | nsolve_grid_size = 16 |
| std::unique_ptr< MLLinOpT< MF > > | ns_linop |
| std::unique_ptr< MLMGT< MF > > | ns_mlmg |
| std::unique_ptr< MF > | ns_sol |
| std::unique_ptr< MF > | ns_rhs |
| std::string | print_ident |
| bool | do_no_sync_gpu = false |
| Vector< MF > | sol |
| Hypre. | |
| Vector< MF > | rhs |
| Vector< int > | sol_is_alias |
| Vector< Vector< MF > > | res |
| First Vector: Amr levels. 0 is the coarest level Second Vector: MG levels. 0 is the finest level. | |
| Vector< Vector< MF > > | cor |
| = rhs - L(sol) | |
| Vector< Vector< MF > > | cor_hold |
| Vector< Vector< MF > > | rescor |
| Vector< double > | timer |
| RT | m_rhsnorm0 = RT(-1.0) |
| RT | m_init_resnorm0 = RT(-1.0) |
| RT | m_final_resnorm0 = RT(-1.0) |
| Vector< int > | m_niters_cg |
| Vector< RT > | m_iter_fine_resnorm0 |
Friends | |
| template<typename T > | |
| class | MLCGSolverT |
| template<typename M > | |
| class | GMRESMLMGT |
| using amrex::MLMGT< MF >::BCMode = typename MLLinOpT<MF>::BCMode |
| using amrex::MLMGT< MF >::BottomSolver = amrex::BottomSolver |
| using amrex::MLMGT< MF >::FAB = typename MLLinOpT<MF>::FAB |
| using amrex::MLMGT< MF >::Location = typename MLLinOpT<MF>::Location |
| using amrex::MLMGT< MF >::MFType = MF |
| using amrex::MLMGT< MF >::RT = typename MLLinOpT<MF>::RT |
|
strong |
|
private |
| amrex::MLMGT< MF >::MLMGT | ( | MLLinOpT< MF > & | a_lp | ) |
|
default |
|
delete |
|
delete |
| void amrex::MLMGT< MF >::actualBottomSolve | ( | ) |
| void amrex::MLMGT< MF >::addInterpCorrection | ( | int | alev, |
| int | mglev | ||
| ) |
| void amrex::MLMGT< MF >::apply | ( | const Vector< MF * > & | out, |
| const Vector< MF * > & | in | ||
| ) |
out = L(in). Note that, if no actual solve is needed, one could turn off multigrid coarsening by constructing a MLLinOp object with an appropriate LPInfo object (e.g., with LPInfo().setMaxCoarseningLevel(0)).
| void amrex::MLMGT< MF >::applyPrecond | ( | const Vector< MF * > & | out, |
| const Vector< MF * > & | in | ||
| ) |
out = L(in) as a preconditioner
| void amrex::MLMGT< MF >::bottomSolve | ( | ) |
| int amrex::MLMGT< MF >::bottomSolveWithCG | ( | MF & | x, |
| const MF & | b, | ||
| typename MLCGSolverT< MF >::Type | type | ||
| ) |
|
private |
| void amrex::MLMGT< MF >::compResidual | ( | const Vector< MF * > & | a_res, |
| const Vector< MF * > & | a_sol, | ||
| const Vector< MF const * > & | a_rhs | ||
| ) |
| void amrex::MLMGT< MF >::computeMLResidual | ( | int | amrlevmax | ) |
| void amrex::MLMGT< MF >::computeResidual | ( | int | alev | ) |
| void amrex::MLMGT< MF >::computeResOfCorrection | ( | int | amrlev, |
| int | mglev | ||
| ) |
| void amrex::MLMGT< MF >::computeResWithCrseCorFineCor | ( | int | falev | ) |
| void amrex::MLMGT< MF >::computeResWithCrseSolFineCor | ( | int | calev, |
| int | falev | ||
| ) |
| void amrex::MLMGT< MF >::decPrintIdentation | ( | ) |
|
inlinenoexcept |
|
inlinenoexcept |
|
inline |
| void amrex::MLMGT< MF >::getEBFluxes | ( | const Vector< MF * > & | a_eb_flux | ) |
| void amrex::MLMGT< MF >::getEBFluxes | ( | const Vector< MF * > & | a_eb_flux, |
| const Vector< MF * > & | a_sol | ||
| ) |
|
inlinenoexcept |
| void amrex::MLMGT< MF >::getFluxes | ( | const Vector< AMF * > & | a_flux, |
| const Vector< AMF * > & | a_sol, | ||
| Location | a_loc = Location::CellCenter |
||
| ) |
| void amrex::MLMGT< MF >::getFluxes | ( | const Vector< AMF * > & | a_flux, |
| Location | a_loc = Location::CellCenter |
||
| ) |
| void amrex::MLMGT< MF >::getFluxes | ( | const Vector< Array< AMF *, 3 > > & | a_flux, |
| const Vector< AMF * > & | a_sol, | ||
| Location | a_loc = Location::FaceCenter |
||
| ) |
| void amrex::MLMGT< MF >::getFluxes | ( | const Vector< Array< AMF *, 3 > > & | a_flux, |
| Location | a_loc = Location::FaceCenter |
||
| ) |
For (alpha * a - beta * (del dot b grad)) phi = rhs, flux means -b grad phi
| void amrex::MLMGT< MF >::getFluxes | ( | std::initializer_list< AMF * > | a_flux, |
| Location | a_loc = Location::CellCenter |
||
| ) |
| void amrex::MLMGT< MF >::getFluxes | ( | std::initializer_list< AMF * > | a_flux, |
| std::initializer_list< AMF * > | a_sol, | ||
| Location | a_loc = Location::CellCenter |
||
| ) |
| void amrex::MLMGT< MF >::getFluxes | ( | std::initializer_list< Array< AMF *, 3 > > | a_flux, |
| Location | a_loc = Location::FaceCenter |
||
| ) |
| void amrex::MLMGT< MF >::getFluxes | ( | std::initializer_list< Array< AMF *, 3 > > | a_flux, |
| std::initializer_list< AMF * > | a_sol, | ||
| Location | a_loc = Location::FaceCenter |
||
| ) |
| void amrex::MLMGT< MF >::getGradSolution | ( | const Vector< Array< AMF *, 3 > > & | a_grad_sol, |
| Location | a_loc = Location::FaceCenter |
||
| ) |
| void amrex::MLMGT< MF >::getGradSolution | ( | std::initializer_list< Array< AMF *, 3 > > | a_grad_sol, |
| Location | a_loc = Location::FaceCenter |
||
| ) |
|
inlinenoexcept |
|
inlinenoexcept |
|
inline |
|
inlinenoexcept |
|
inlinenoexcept |
|
inlinenoexcept |
|
inline |
| void amrex::MLMGT< MF >::incPrintIdentation | ( | ) |
| void amrex::MLMGT< MF >::interpCorrection | ( | int | alev | ) |
| void amrex::MLMGT< MF >::interpCorrection | ( | int | alev, |
| int | mglev | ||
| ) |
| void amrex::MLMGT< MF >::makeSolvable | ( | ) |
| void amrex::MLMGT< MF >::makeSolvable | ( | int | amrlev, |
| int | mglev, | ||
| MF & | mf | ||
| ) |
| void amrex::MLMGT< MF >::mgFcycle | ( | ) |
| void amrex::MLMGT< MF >::mgVcycle | ( | int | amrlev, |
| int | mglev | ||
| ) |
| void amrex::MLMGT< MF >::miniCycle | ( | int | amrlev | ) |
| auto amrex::MLMGT< MF >::MLResNormInf | ( | int | alevmax, |
| bool | local = false |
||
| ) |
| auto amrex::MLMGT< MF >::MLRhsNormInf | ( | bool | local = false | ) |
| void amrex::MLMGT< MF >::NSolve | ( | MLMGT< MF > & | a_solver, |
| MF & | a_sol, | ||
| MF & | a_rhs | ||
| ) |
|
inlinenoexcept |
| void amrex::MLMGT< MF >::oneIter | ( | int | iter | ) |
|
delete |
|
delete |
| auto amrex::MLMGT< MF >::precond | ( | Vector< MF * > const & | a_sol, |
| Vector< MF const * > const & | a_rhs, | ||
| RT | a_tol_rel, | ||
| RT | a_tol_abs | ||
| ) |
| void amrex::MLMGT< MF >::prepareForFluxes | ( | Vector< MF const * > const & | a_sol | ) |
| void amrex::MLMGT< MF >::prepareForNSolve | ( | ) |
| void amrex::MLMGT< MF >::prepareForSolve | ( | Vector< AMF * > const & | a_sol, |
| Vector< AMF const * > const & | a_rhs | ||
| ) |
| void amrex::MLMGT< MF >::prepareLinOp | ( | ) |
| void amrex::MLMGT< MF >::preparePrecond | ( | ) |
| auto amrex::MLMGT< MF >::ResNormInf | ( | int | alev, |
| bool | local = false |
||
| ) |
|
noexcept |
|
inlinenoexcept |
|
inlinenoexcept |
|
inlinenoexcept |
|
inlinenoexcept |
|
inlinenoexcept |
|
inlinenoexcept |
|
inlinenoexcept |
|
inlinenoexcept |
|
inlinenoexcept |
|
inlinenoexcept |
|
inlinenoexcept |
|
inlinenoexcept |
|
inlinenoexcept |
|
inlinenoexcept |
|
inlinenoexcept |
|
inlinenoexcept |
|
inlinenoexcept |
|
inlinenoexcept |
|
inlinenoexcept |
|
inlinenoexcept |
|
inlinenoexcept |
| auto amrex::MLMGT< MF >::solve | ( | const Vector< AMF * > & | a_sol, |
| const Vector< AMF const * > & | a_rhs, | ||
| RT | a_tol_rel, | ||
| RT | a_tol_abs, | ||
| const char * | checkpoint_file | ||
| ) | -> RT |
| RT amrex::MLMGT< MF >::solve | ( | const Vector< AMF * > & | a_sol, |
| const Vector< AMF const * > & | a_rhs, | ||
| RT | a_tol_rel, | ||
| RT | a_tol_abs, | ||
| const char * | checkpoint_file = nullptr |
||
| ) |
| auto amrex::MLMGT< MF >::solve | ( | std::initializer_list< AMF * > | a_sol, |
| std::initializer_list< AMF const * > | a_rhs, | ||
| RT | a_tol_rel, | ||
| RT | a_tol_abs, | ||
| const char * | checkpoint_file | ||
| ) | -> RT |
| RT amrex::MLMGT< MF >::solve | ( | std::initializer_list< AMF * > | a_sol, |
| std::initializer_list< AMF const * > | a_rhs, | ||
| RT | a_tol_rel, | ||
| RT | a_tol_abs, | ||
| const char * | checkpoint_file = nullptr |
||
| ) |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
= rhs - L(sol)
L(cor) = res
|
private |
|
private |
|
private |
|
private |
N Solve.
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
pre
|
private |
post
|
private |
additional smoothing after bottom cg solver
|
private |
when smoother is used as bottom solver
|
private |
|
private |
|
private |
|
private |
= res - L(cor) Residual of the correction form
|
private |
Copy of original rhs L(sol) = rhs
|
private |
PETSc
To avoid confusion, terms like sol, cor, rhs, res, ... etc. are in the frame of the original equation, not the correction form Might be alias to argument a_sol
|
private |
|
private |
|
private |
|
private |
|
private |