Multigrid solver for linear systems. More...
#include <multigrid.hpp>
Public Types | |
| using | OperatorType = OperatorT |
| Operator type to be solved. | |
| using | ProlongationType = ProlongationT |
| Prolongation operator type. | |
| using | RestrictionType = RestrictionT |
| Restriction operator type. | |
| using | SmootherType = SmootherT |
| Smoother type. | |
| using | CoarseGridSolverType = CoarseGridSolverT |
| Coarse grid solver type. | |
| using | SolutionVectorType = SrcOf< OperatorType > |
| Solution vector type. | |
| using | RHSVectorType = DstOf< OperatorType > |
| Right-hand side vector type. | |
| using | ScalarType = SolutionVectorType::ScalarType |
| Scalar type for computations. | |
Public Member Functions | |
| Multigrid (const std::vector< ProlongationType > &P_additive, const std::vector< RestrictionType > &R, const std::vector< OperatorT > &A_c, const std::vector< SolutionVectorType > &tmp_r, const std::vector< SolutionVectorType > &tmp_e, const std::vector< SolutionVectorType > &tmp, const std::vector< SmootherType > &smoothers_pre, const std::vector< SmootherType > &smoothers_post, const CoarseGridSolverType &coarse_grid_solver, int num_cycles, ScalarType relative_residual_threshold) | |
| Construct a multigrid solver. | |
| Multigrid (const std::vector< ProlongationType > &P_additive, const std::vector< RestrictionType > &R, const std::vector< OperatorT > &A_c, const std::vector< SolutionVectorType > &tmp_r, const std::vector< SolutionVectorType > &tmp_e, const std::vector< SolutionVectorType > &tmp, const std::vector< SmootherType > &smoothers_pre, const std::vector< SmootherType > &smoothers_post, const CoarseGridSolverType &coarse_grid_solver, int num_cycles, ScalarType relative_residual_threshold, std::vector< RedistributeT > redistribute_down, std::vector< SolutionVectorType > tmp_r_fine, std::vector< SolutionVectorType > tmp_e_fine) | |
| Construct an agglomerated multigrid solver. | |
| void | set_tag (const std::string &tag) |
| Set a tag string for statistics output. | |
| void | collect_statistics (const std::shared_ptr< util::Table > &statistics) |
| Collect statistics in a shared table. | |
| void | solve_impl (OperatorType &A, SolutionVectorType &x, const RHSVectorType &b) |
| Solve the linear system using multigrid cycles. Calls the recursive V-cycle and updates statistics. | |
Multigrid solver for linear systems.
Satisfies the SolverLike concept (see solver.hpp). Supports arbitrary operators, prolongation/restriction, smoothers, and coarse grid solvers. Implements recursive V-cycle multigrid.
Optionally supports hierarchical agglomeration: per level, a Redistribute operator moves the restricted residual from the fine level's communicator to a (smaller) coarse level communicator, and the coarse correction back up. Inactive coarse ranks skip the recursive descent but still participate in the collectives inside Redistribute::apply / apply_transpose.
| OperatorT | Operator type (must satisfy OperatorLike). |
| ProlongationT | Prolongation operator type (must satisfy OperatorLike). |
| RestrictionT | Restriction operator type (must satisfy OperatorLike). |
| SmootherT | Smoother type (must satisfy SolverLike). |
| CoarseGridSolverT | Coarse grid solver type (must satisfy SolverLike). |
| RedistributeT | Optional; defaults to NoRedistribute. Pass the concrete Redistribute<GridDataType> type when using agglomeration. |
| using terra::linalg::solvers::Multigrid< OperatorT, ProlongationT, RestrictionT, SmootherT, CoarseGridSolverT, RedistributeT >::CoarseGridSolverType = CoarseGridSolverT |
Coarse grid solver type.
| using terra::linalg::solvers::Multigrid< OperatorT, ProlongationT, RestrictionT, SmootherT, CoarseGridSolverT, RedistributeT >::OperatorType = OperatorT |
Operator type to be solved.
| using terra::linalg::solvers::Multigrid< OperatorT, ProlongationT, RestrictionT, SmootherT, CoarseGridSolverT, RedistributeT >::ProlongationType = ProlongationT |
Prolongation operator type.
| using terra::linalg::solvers::Multigrid< OperatorT, ProlongationT, RestrictionT, SmootherT, CoarseGridSolverT, RedistributeT >::RestrictionType = RestrictionT |
Restriction operator type.
| using terra::linalg::solvers::Multigrid< OperatorT, ProlongationT, RestrictionT, SmootherT, CoarseGridSolverT, RedistributeT >::RHSVectorType = DstOf< OperatorType > |
Right-hand side vector type.
| using terra::linalg::solvers::Multigrid< OperatorT, ProlongationT, RestrictionT, SmootherT, CoarseGridSolverT, RedistributeT >::ScalarType = SolutionVectorType::ScalarType |
Scalar type for computations.
| using terra::linalg::solvers::Multigrid< OperatorT, ProlongationT, RestrictionT, SmootherT, CoarseGridSolverT, RedistributeT >::SmootherType = SmootherT |
Smoother type.
| using terra::linalg::solvers::Multigrid< OperatorT, ProlongationT, RestrictionT, SmootherT, CoarseGridSolverT, RedistributeT >::SolutionVectorType = SrcOf< OperatorType > |
Solution vector type.
|
inline |
Construct a multigrid solver.
Vector ordering of arguments always goes from the coarsest level (index 0) to the finest.
| P_additive | Prolongation operators for each coarse level. Size must match the number of levels - 1. Must be additive prolongation operators, i.e., apply( P, x, y )
void apply(LinearForm &L, typename LinearForm::DstVectorType &dst) Apply a linear form and write to a destination vector. Definition linear_form.hpp:37 |
| R | Restriction operators for each coarse level. |
| A_c | Coarse grid operators for each coarse level. |
| tmp_r | Temporary residual vectors for each coarse level. |
| tmp_e | Temporary error vectors for each coarse level. |
| tmp | Temporary workspace vectors for each level (including the finest level). |
| smoothers_pre | Pre-smoothers for each level (including the finest level). |
| smoothers_post | Post-smoothers for each level (including the finest level). |
| coarse_grid_solver | Coarse grid solver. |
| num_cycles | Number of multigrid cycles to perform. |
| relative_residual_threshold | Relative residual threshold for stopping. |
|
inline |
Construct an agglomerated multigrid solver.
Same as the non-agglomerated constructor, plus:
redistribute_down: Redistribute operators, one per coarse descent step. redistribute_down[L-1] moves a vector from level L's comm to level L-1's comm. apply_transpose is used on the way back up.tmp_r_fine, tmp_e_fine: per-coarse-level scratch vectors allocated on the fine level's communicator. Restriction output goes into tmp_r_fine; prolongation input is read from tmp_e_fine. tmp_r, tmp_e must carry the coarse level's sub-comm; on ranks not part of the sub-comm they are empty (zero-subdomain) vectors.
|
inline |
Collect statistics in a shared table.
| statistics | Shared pointer to statistics table. |
|
inline |
Set a tag string for statistics output.
| tag | Tag string. |
|
inline |
Solve the linear system using multigrid cycles. Calls the recursive V-cycle and updates statistics.
| A | Operator (matrix). |
| x | Solution vector (output). |
| b | Right-hand side vector (input). |