AlgebraicMultigrid 0.1
C++ algebraic multigrid.
|
Interface for multigrid solver (via solve
) that iteratively performs standard vcycle
s.
More...
#include <multigrid.hpp>
Public Member Functions | |
~Multigrid () | |
Multigrid ()=delete | |
Multigrid (AMG::InterpolatorBase< EleType > *interpolator_, AMG::SmootherBase< EleType > *smoother_, const Eigen::SparseMatrix< EleType > &A, const Eigen::Matrix< EleType, -1, 1 > &b, size_t n_levels_, EleType tolerance_=1e-9, size_t compute_error_every_n_iters_=10, size_t n_iters_=100) | |
Construct a new Multigrid object. | |
void | vcycle () |
A single multigrid vcycle. | |
const Eigen::Matrix< EleType, -1, 1 > & | solve () |
Solve a linear system of equations via vcycle. | |
const Eigen::SparseMatrix< EleType > & | get_coefficient_matrix (size_t level) const |
const Eigen::Matrix< EleType, -1, 1 > & | get_soln (size_t level) const |
const Eigen::Matrix< EleType, -1, 1 > & | get_rhs (size_t level) const |
const size_t | get_n_dofs (size_t level) const |
const EleType | get_tolerance () const |
void | display_error_on () |
void | display_error_off () |
Private Member Functions | |
size_t | n_H_dofs_from_n_h_dofs (size_t h_dofs) |
Return number of next coarsest dofs given finer dofs. | |
Private Attributes | |
InterpolatorBase< EleType > * | interpolator |
Interpolator providing restriction and prolongation operations. | |
AMG::SmootherBase< EleType > * | smoother |
Eigen::SimplicialLDLT< Eigen::SparseMatrix< EleType > > | coarse_direct_solver |
EleType | tolerance |
Tolerance below which a smoother is considered to have converged. | |
size_t | compute_error_every_n_iters |
Compute the error every n iterations during smoothing. | |
size_t | n_iters |
Maximum number of iterations before smoothing termination. | |
size_t | n_fine_nodes |
size_t | n_levels |
Number of multigrid levels where level 0 is finest. | |
size_t | finest_grid_ix = 0 |
Index of the finest grid. | |
size_t | coarsest_grid_ix |
Index of the coarsest grid. | |
std::vector< size_t > | level_to_n_dofs |
Map multigrid level to number of degrees of freedom for that level. | |
std::vector< Eigen::SparseMatrix< EleType > > | level_to_coefficient_matrix |
Map multigrid level to coefficient matrix A. | |
std::vector< Eigen::Matrix< EleType, -1, 1 > > | level_to_rhs |
Map multigrid level to right hand side (forcing) vector. | |
std::vector< Eigen::Matrix< EleType, -1, 1 > > | level_to_soln |
Map multigrid level to solution (u) vector. | |
std::vector< Eigen::Matrix< EleType, -1, 1 > > | level_to_residual |
Map multigrid level to residual (e) vector. | |
bool | display_error {false} |
Interface for multigrid solver (via solve
) that iteratively performs standard vcycle
s.
EleType |
|
inline |
|
delete |
|
inline |
Construct a new Multigrid object.
interpolator_ | Provide prolongation and restriction operators. |
smoother_ | Smooth high frequency errors in post/pre relaxation. |
A | Coefficients matrix for discretized governing equations. |
b | Right hand side of linear system Au = b . |
n_levels_ | Desired number of levels where level 0 is finest level. |
tolerance_ | Minimum tolerance residual must be leq to converge. |
compute_error_every_n_iters_ | Compute convergence on this interval. |
n_iters_ | Max number of iterations for stopping solution via vcycles. |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inlineprivate |
Return number of next coarsest dofs given finer dofs.
The next coarsest dofs formula follows from the following (see ref [1]):
References:
[1] : Briggs2000. "A Multigrid Tutorial, 2ed". pp. 34.
|
inline |
Solve a linear system of equations via vcycle.
|
inline |
A single multigrid vcycle.
Ref [1] was used for the main algorithm (non-recursive) while refs [2-4] were used to verify that coarse solutions need to be zeroed at each vcycle.
References:
[1] : Algebraic Multigrid from amgcl
[2] : Kostler2008. "Multigrid HowTo: simple Multigrid solver in C++ in less than 200 lines of code"
[3] : Pawar2019. "CFD Julia: A Learning Module Structuring an Introductory Course on Computational Fluid Dynamics "
|
private |
|
private |
Index of the coarsest grid.
|
private |
Compute the error every n
iterations during smoothing.
|
private |
|
private |
Index of the finest grid.
|
private |
Interpolator providing restriction and prolongation operations.
|
private |
Map multigrid level to coefficient matrix A.
|
private |
Map multigrid level to number of degrees of freedom for that level.
Note that n_dofs==n_nodes*n_nodes
|
private |
Map multigrid level to residual (e) vector.
|
private |
Map multigrid level to right hand side (forcing) vector.
|
private |
Map multigrid level to solution (u) vector.
|
private |
|
private |
Maximum number of iterations before smoothing termination.
|
private |
Number of multigrid levels where level 0 is finest.
|
private |
|
private |
Tolerance below which a smoother is considered to have converged.