Dune Core Modules (2.10.0)
An Algebraic Multigrid based on Agglomeration that saves memory bandwidth. More...
Classes | |
class | Dune::Amg::FastAMG< M, X, PI, A > |
A fast (sequential) algebraic multigrid based on agglomeration that saves memory bandwidth. More... | |
Typedefs | |
typedef M | Dune::Amg::FastAMG< M, X, PI, A >::Operator |
The matrix operator type. | |
typedef PI | Dune::Amg::FastAMG< M, X, PI, A >::ParallelInformation |
The type of the parallel information. Either OwnerOverlapCommunication or another type describing the parallel data distribution and providing communication methods. | |
typedef MatrixHierarchy< M, ParallelInformation, A > | Dune::Amg::FastAMG< M, X, PI, A >::OperatorHierarchy |
The operator hierarchy type. | |
typedef OperatorHierarchy::ParallelInformationHierarchy | Dune::Amg::FastAMG< M, X, PI, A >::ParallelInformationHierarchy |
The parallal data distribution hierarchy type. | |
typedef X | Dune::Amg::FastAMG< M, X, PI, A >::Domain |
The domain type. | |
typedef X | Dune::Amg::FastAMG< M, X, PI, A >::Range |
The range type. | |
typedef InverseOperator< X, X > | Dune::Amg::FastAMG< M, X, PI, A >::CoarseSolver |
the type of the coarse solver. | |
Functions | |
Dune::Amg::FastAMG< M, X, PI, A >::FastAMG (OperatorHierarchy &matrices, CoarseSolver &coarseSolver, const Parameters &parms, bool symmetric=true) | |
Construct a new amg with a specific coarse solver. More... | |
template<class C > | |
Dune::Amg::FastAMG< M, X, PI, A >::FastAMG (std::shared_ptr< const Operator > fineOperator, const C &criterion, const Parameters &parms=Parameters(), bool symmetric=true, const ParallelInformation &pinfo=ParallelInformation()) | |
Construct an AMG with an inexact coarse solver based on the smoother. More... | |
template<class C > | |
Dune::Amg::FastAMG< M, X, PI, A >::FastAMG (const Operator &fineOperator, const C &criterion, const Parameters &parms=Parameters(), bool symmetric=true, const ParallelInformation &pinfo=ParallelInformation()) | |
Construct an AMG with an inexact coarse solver based on the smoother. More... | |
Dune::Amg::FastAMG< M, X, PI, A >::FastAMG (const FastAMG &amg) | |
Copy constructor. | |
void | Dune::Amg::FastAMG< M, X, PI, A >::pre (Domain &x, Range &b) |
Prepare the preconditioner. More... | |
void | Dune::Amg::FastAMG< M, X, PI, A >::apply (Domain &v, const Range &d) |
Apply one step of the preconditioner to the system A(v)=d. More... | |
virtual SolverCategory::Category | Dune::Amg::FastAMG< M, X, PI, A >::category () const |
Category of the preconditioner (see SolverCategory::Category) | |
void | Dune::Amg::FastAMG< M, X, PI, A >::post (Domain &x) |
Clean up. More... | |
template<class A1 > | |
void | Dune::Amg::FastAMG< M, X, PI, A >::getCoarsestAggregateNumbers (std::vector< std::size_t, A1 > &cont) |
Get the aggregate number of each unknown on the coarsest level. More... | |
void | Dune::Amg::FastAMG< M, X, PI, A >::recalculateHierarchy () |
Recalculate the matrix hierarchy. More... | |
bool | Dune::Amg::FastAMG< M, X, PI, A >::usesDirectCoarseLevelSolver () const |
Check whether the coarse solver used is a direct solver. More... | |
Variables | |
OperatorHierarchy::ParallelMatrixHierarchy::ConstIterator | Dune::Amg::FastAMG< M, X, PI, A >::LevelContext::matrix |
The iterator over the matrices. | |
ParallelInformationHierarchy::Iterator | Dune::Amg::FastAMG< M, X, PI, A >::LevelContext::pinfo |
The iterator over the parallel information. | |
OperatorHierarchy::RedistributeInfoList::const_iterator | Dune::Amg::FastAMG< M, X, PI, A >::LevelContext::redist |
The iterator over the redistribution information. | |
OperatorHierarchy::AggregatesMapList::const_iterator | Dune::Amg::FastAMG< M, X, PI, A >::LevelContext::aggregates |
The iterator over the aggregates maps. | |
Hierarchy< Domain, A >::Iterator | Dune::Amg::FastAMG< M, X, PI, A >::LevelContext::lhs |
The iterator over the left hand side. | |
Hierarchy< Domain, A >::Iterator | Dune::Amg::FastAMG< M, X, PI, A >::LevelContext::residual |
The iterator over the residuals. | |
Hierarchy< Range, A >::Iterator | Dune::Amg::FastAMG< M, X, PI, A >::LevelContext::rhs |
The iterator over the right hand sided. | |
std::size_t | Dune::Amg::FastAMG< M, X, PI, A >::LevelContext::level |
The level index. | |
Detailed Description
An Algebraic Multigrid based on Agglomeration that saves memory bandwidth.
Function Documentation
◆ apply()
void Dune::Amg::FastAMG< M, X, PI, A >::apply | ( | Domain & | v, |
const Range & | d | ||
) |
Apply one step of the preconditioner to the system A(v)=d.
On entry v=0 and d=b-A(x) (although this might not be computed in that way. On exit v contains the update, i.e one step computes \( v = M^{-1} d \) where \( M \) is the approximate inverse of the operator \( A \) characterizing the preconditioner.
- Parameters
-
[out] v The update to be computed d The current defect.
◆ FastAMG() [1/3]
|
inline |
Construct an AMG with an inexact coarse solver based on the smoother.
As coarse solver a preconditioned CG method with the smoother as preconditioner will be used. The matrix hierarchy is built automatically.
- Parameters
-
fineOperator The operator on the fine level. criterion The criterion describing the coarsening strategy. E. g. SymmetricCriterion or UnsymmetricCriterion, and providing the parameters. parms The parameters for the AMG.
NOTE: This constructor uses the default constructed ParallelInformation
.
◆ FastAMG() [2/3]
Dune::Amg::FastAMG< M, X, PI, A >::FastAMG | ( | OperatorHierarchy & | matrices, |
CoarseSolver & | coarseSolver, | ||
const Parameters & | parms, | ||
bool | symmetric = true |
||
) |
Construct a new amg with a specific coarse solver.
- Parameters
-
matrices The already set-up matrix hierarchy. coarseSolver The set up solver to use on the coarse grid, must match the coarse matrix in the matrix hierarchy. parms The parameters for the AMG.
◆ FastAMG() [3/3]
Dune::Amg::FastAMG< M, X, PI, A >::FastAMG | ( | std::shared_ptr< const Operator > | fineOperator, |
const C & | criterion, | ||
const Parameters & | parms = Parameters() , |
||
bool | symmetric = true , |
||
const ParallelInformation & | pinfo = ParallelInformation() |
||
) |
Construct an AMG with an inexact coarse solver based on the smoother.
As coarse solver a preconditioned CG method with the smoother as preconditioner will be used. The matrix hierarchy is built automatically.
- Parameters
-
fineOperator The operator on the fine level. criterion The criterion describing the coarsening strategy. E. g. SymmetricCriterion or UnsymmetricCriterion, and providing the parameters. parms The parameters for the AMG. pinfo The information about the parallel distribution of the data.
◆ getCoarsestAggregateNumbers()
void Dune::Amg::FastAMG< M, X, PI, A >::getCoarsestAggregateNumbers | ( | std::vector< std::size_t, A1 > & | cont | ) |
Get the aggregate number of each unknown on the coarsest level.
- Parameters
-
cont The random access container to store the numbers in.
◆ post()
|
virtual |
Clean up.
This method is called after the last apply call for the linear system to be solved. Memory may be deallocated safely here. x is the solution of the linear equation.
- Parameters
-
x The right hand side of the equation.
Implements Dune::Preconditioner< X, X >.
◆ pre()
void Dune::Amg::FastAMG< M, X, PI, A >::pre | ( | Domain & | x, |
Range & | b | ||
) |
Prepare the preconditioner.
A solver solves a linear operator equation A(x)=b by applying one or several steps of the preconditioner. The method pre() is called before the first apply operation. b and x are right hand side and solution vector of the linear system respectively. It may. e.g., scale the system, allocate memory or compute a (I)LU decomposition. Note: The ILU decomposition could also be computed in the constructor or with a separate method of the derived method if several linear systems with the same matrix are to be solved.
- Note
- if a preconditioner is copied (e.g. for a second thread) again the pre() method has to be called to ensure proper memory management.
- Parameters
-
x The left hand side of the equation. b The right hand side of the equation.
◆ recalculateHierarchy()
|
inline |
Recalculate the matrix hierarchy.
It is assumed that the coarsening for the changed fine level matrix would yield the same aggregates. In this case it suffices to recalculate all the Galerkin products for the matrices of the coarser levels.
◆ usesDirectCoarseLevelSolver()
bool Dune::Amg::FastAMG< M, X, PI, A >::usesDirectCoarseLevelSolver |
Check whether the coarse solver used is a direct solver.
- Returns
- True if the coarse level solver is a direct solver.