AlgebraicMultigrid 0.1
C++ algebraic multigrid.
Loading...
Searching...
No Matches
AMG::Multigrid< EleType > Class Template Reference

Interface for multigrid solver (via solve) that iteratively performs standard vcycles. 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}
 

Detailed Description

template<class EleType>
class AMG::Multigrid< EleType >

Interface for multigrid solver (via solve) that iteratively performs standard vcycles.

Template Parameters
EleType

Constructor & Destructor Documentation

◆ ~Multigrid()

template<class EleType >
AMG::Multigrid< EleType >::~Multigrid ( )
inline

◆ Multigrid() [1/2]

template<class EleType >
AMG::Multigrid< EleType >::Multigrid ( )
delete

◆ Multigrid() [2/2]

template<class EleType >
AMG::Multigrid< EleType >::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 )
inline

Construct a new Multigrid object.

Parameters
interpolator_Provide prolongation and restriction operators.
smoother_Smooth high frequency errors in post/pre relaxation.
ACoefficients matrix for discretized governing equations.
bRight 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.

Member Function Documentation

◆ display_error_off()

template<class EleType >
void AMG::Multigrid< EleType >::display_error_off ( )
inline

◆ display_error_on()

template<class EleType >
void AMG::Multigrid< EleType >::display_error_on ( )
inline

◆ get_coefficient_matrix()

template<class EleType >
const Eigen::SparseMatrix< EleType > & AMG::Multigrid< EleType >::get_coefficient_matrix ( size_t level) const
inline

◆ get_n_dofs()

template<class EleType >
const size_t AMG::Multigrid< EleType >::get_n_dofs ( size_t level) const
inline

◆ get_rhs()

template<class EleType >
const Eigen::Matrix< EleType, -1, 1 > & AMG::Multigrid< EleType >::get_rhs ( size_t level) const
inline

◆ get_soln()

template<class EleType >
const Eigen::Matrix< EleType, -1, 1 > & AMG::Multigrid< EleType >::get_soln ( size_t level) const
inline

◆ get_tolerance()

template<class EleType >
const EleType AMG::Multigrid< EleType >::get_tolerance ( ) const
inline

◆ n_H_dofs_from_n_h_dofs()

template<class EleType >
size_t AMG::Multigrid< EleType >::n_H_dofs_from_n_h_dofs ( size_t h_dofs)
inlineprivate

Return number of next coarsest dofs given finer dofs.

The next coarsest dofs formula follows from the following (see ref [1]):

nc = n/2 - 1
nf = n - 1
# nc in terms of nf
nc = (nf + 1)/2 - 1

References:

[1] : Briggs2000. "A Multigrid Tutorial, 2ed". pp. 34.

Returns
size_t

◆ solve()

template<class EleType >
const Eigen::Matrix< EleType, -1, 1 > & AMG::Multigrid< EleType >::solve ( )
inline

Solve a linear system of equations via vcycle.

◆ vcycle()

template<class EleType >
void AMG::Multigrid< EleType >::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 "

[4] : AlgebraicMultigrid.jl solve

Member Data Documentation

◆ coarse_direct_solver

template<class EleType >
Eigen::SimplicialLDLT<Eigen::SparseMatrix<EleType> > AMG::Multigrid< EleType >::coarse_direct_solver
private

◆ coarsest_grid_ix

template<class EleType >
size_t AMG::Multigrid< EleType >::coarsest_grid_ix
private

Index of the coarsest grid.

◆ compute_error_every_n_iters

template<class EleType >
size_t AMG::Multigrid< EleType >::compute_error_every_n_iters
private

Compute the error every n iterations during smoothing.

◆ display_error

template<class EleType >
bool AMG::Multigrid< EleType >::display_error {false}
private

◆ finest_grid_ix

template<class EleType >
size_t AMG::Multigrid< EleType >::finest_grid_ix = 0
private

Index of the finest grid.

◆ interpolator

template<class EleType >
InterpolatorBase<EleType>* AMG::Multigrid< EleType >::interpolator
private

Interpolator providing restriction and prolongation operations.

◆ level_to_coefficient_matrix

template<class EleType >
std::vector<Eigen::SparseMatrix<EleType> > AMG::Multigrid< EleType >::level_to_coefficient_matrix
private

Map multigrid level to coefficient matrix A.

◆ level_to_n_dofs

template<class EleType >
std::vector<size_t> AMG::Multigrid< EleType >::level_to_n_dofs
private

Map multigrid level to number of degrees of freedom for that level.

Note that n_dofs==n_nodes*n_nodes

◆ level_to_residual

template<class EleType >
std::vector<Eigen::Matrix<EleType, -1, 1> > AMG::Multigrid< EleType >::level_to_residual
private

Map multigrid level to residual (e) vector.

◆ level_to_rhs

template<class EleType >
std::vector<Eigen::Matrix<EleType, -1, 1> > AMG::Multigrid< EleType >::level_to_rhs
private

Map multigrid level to right hand side (forcing) vector.

◆ level_to_soln

template<class EleType >
std::vector<Eigen::Matrix<EleType, -1, 1> > AMG::Multigrid< EleType >::level_to_soln
private

Map multigrid level to solution (u) vector.

◆ n_fine_nodes

template<class EleType >
size_t AMG::Multigrid< EleType >::n_fine_nodes
private

◆ n_iters

template<class EleType >
size_t AMG::Multigrid< EleType >::n_iters
private

Maximum number of iterations before smoothing termination.

◆ n_levels

template<class EleType >
size_t AMG::Multigrid< EleType >::n_levels
private

Number of multigrid levels where level 0 is finest.

◆ smoother

template<class EleType >
AMG::SmootherBase<EleType>* AMG::Multigrid< EleType >::smoother
private

◆ tolerance

template<class EleType >
EleType AMG::Multigrid< EleType >::tolerance
private

Tolerance below which a smoother is considered to have converged.


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