dune-istl  2.4.1
Classes | Typedefs | Enumerations | Enumerator | Functions | Variables

An Algebraic Multigrid based on Agglomeration that saves memory bandwidth. More...

Collaboration diagram for Fast (sequential) Algebraic Multigrid:

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. More...
 
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. More...
 
typedef MatrixHierarchy< M,
ParallelInformation, A
Dune::Amg::FastAMG< M, X, PI, A >::OperatorHierarchy
 The operator hierarchy type. More...
 
typedef
OperatorHierarchy::ParallelInformationHierarchy 
Dune::Amg::FastAMG< M, X, PI, A >::ParallelInformationHierarchy
 The parallal data distribution hierarchy type. More...
 
typedef X Dune::Amg::FastAMG< M, X, PI, A >::Domain
 The domain type. More...
 
typedef X Dune::Amg::FastAMG< M, X, PI, A >::Range
 The range type. More...
 
typedef InverseOperator< X, X > Dune::Amg::FastAMG< M, X, PI, A >::CoarseSolver
 the type of the coarse solver. More...
 

Enumerations

enum  { Dune::Amg::FastAMG< M, X, PI, A >::category = SolverCategory::sequential }
 

Functions

 Dune::Amg::FastAMG< M, X, PI, A >::FastAMG (const 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 (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. More...
 
 Dune::Amg::FastAMG< M, X, PI, A >::~FastAMG ()
 
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...
 
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...
 
std::size_t Dune::Amg::FastAMG< M, X, PI, A >::levels ()
 
std::size_t Dune::Amg::FastAMG< M, X, PI, A >::maxlevels ()
 
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. More...
 
ParallelInformationHierarchy::Iterator Dune::Amg::FastAMG< M, X, PI, A >::LevelContext::pinfo
 The iterator over the parallel information. More...
 
OperatorHierarchy::RedistributeInfoList::const_iterator Dune::Amg::FastAMG< M, X, PI, A >::LevelContext::redist
 The iterator over the redistribution information. More...
 
OperatorHierarchy::AggregatesMapList::const_iterator Dune::Amg::FastAMG< M, X, PI, A >::LevelContext::aggregates
 The iterator over the aggregates maps. More...
 
Hierarchy< Domain, A >::Iterator Dune::Amg::FastAMG< M, X, PI, A >::LevelContext::lhs
 The iterator over the left hand side. More...
 
Hierarchy< Domain, A >::Iterator Dune::Amg::FastAMG< M, X, PI, A >::LevelContext::residual
 The iterator over the residuals. More...
 
Hierarchy< Range, A >::Iterator Dune::Amg::FastAMG< M, X, PI, A >::LevelContext::rhs
 The iterator over the right hand sided. More...
 
std::size_t Dune::Amg::FastAMG< M, X, PI, A >::LevelContext::level
 The level index. More...
 

Detailed Description

An Algebraic Multigrid based on Agglomeration that saves memory bandwidth.

Typedef Documentation

template<class M , class X , class PI = SequentialInformation, class A = std::allocator<X>>
typedef InverseOperator<X,X> Dune::Amg::FastAMG< M, X, PI, A >::CoarseSolver

the type of the coarse solver.

template<class M , class X , class PI = SequentialInformation, class A = std::allocator<X>>
typedef X Dune::Amg::FastAMG< M, X, PI, A >::Domain

The domain type.

template<class M , class X , class PI = SequentialInformation, class A = std::allocator<X>>
typedef M Dune::Amg::FastAMG< M, X, PI, A >::Operator

The matrix operator type.

template<class M , class X , class PI = SequentialInformation, class A = std::allocator<X>>
typedef MatrixHierarchy<M, ParallelInformation, A> Dune::Amg::FastAMG< M, X, PI, A >::OperatorHierarchy

The operator hierarchy type.

template<class M , class X , class PI = SequentialInformation, class A = std::allocator<X>>
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.

template<class M , class X , class PI = SequentialInformation, class A = std::allocator<X>>
typedef OperatorHierarchy::ParallelInformationHierarchy Dune::Amg::FastAMG< M, X, PI, A >::ParallelInformationHierarchy

The parallal data distribution hierarchy type.

template<class M , class X , class PI = SequentialInformation, class A = std::allocator<X>>
typedef X Dune::Amg::FastAMG< M, X, PI, A >::Range

The range type.

Enumeration Type Documentation

template<class M , class X , class PI = SequentialInformation, class A = std::allocator<X>>
anonymous enum
Enumerator
category 

The solver category.

Function Documentation

template<class M , class X , class PI , class A >
void Dune::Amg::FastAMG< M, X, PI, A >::apply ( Domain v,
const Range d 
)
virtual

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]vThe update to be computed
dThe current defect.

Implements Dune::Preconditioner< X, X >.

template<class M , class X , class PI , class A >
Dune::Amg::FastAMG< M, X, PI, A >::FastAMG ( const OperatorHierarchy matrices,
CoarseSolver coarseSolver,
const Parameters parms,
bool  symmetric = true 
)

Construct a new amg with a specific coarse solver.

Parameters
matricesThe already set up matix hierarchy.
coarseSolverThe set up solver to use on the coarse grid, must match the coarse matrix in the matrix hierachy.
parmsThe parameters for the AMG.
template<class M , class X , class PI , class A >
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.

As coarse solver a preconditioned CG method with the smoother as preconditioner will be used. The matrix hierarchy is built automatically.

Parameters
fineOperatorThe operator on the fine level.
criterionThe criterion describing the coarsening strategy. E. g. SymmetricCriterion or UnsymmetricCriterion, and providing the parameters.
parmsThe parameters for the AMG.
pinfoThe information about the parallel distribution of the data.
template<class M , class X , class PI , class A >
Dune::Amg::FastAMG< M, X, PI, A >::FastAMG ( const FastAMG< M, X, PI, A > &  amg)

Copy constructor.

template<class M , class X , class PI , class A >
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.

Parameters
contThe random access container to store the numbers in.
template<class M , class X , class PI , class A >
std::size_t Dune::Amg::FastAMG< M, X, PI, A >::levels ( )
template<class M , class X , class PI , class A >
std::size_t Dune::Amg::FastAMG< M, X, PI, A >::maxlevels ( )
template<class M , class X , class PI , class A >
void Dune::Amg::FastAMG< M, X, PI, A >::post ( Domain x)
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
xThe right hand side of the equation.

Implements Dune::Preconditioner< X, X >.

template<class M , class X , class PI , class A >
void Dune::Amg::FastAMG< M, X, PI, A >::pre ( Domain x,
Range b 
)
virtual

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.

Parameters
xThe left hand side of the equation.
bThe right hand side of the equation.

Implements Dune::Preconditioner< X, X >.

References Dune::Matrix< T, A >::begin(), col, Dune::Matrix< T, A >::end(), Dune::Amg::Hierarchy< T, A >::finest(), mat, and row.

template<class M , class X , class PI = SequentialInformation, class A = std::allocator<X>>
void Dune::Amg::FastAMG< M, X, PI, A >::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.

template<class M , class X , class PI , class A >
bool Dune::Amg::FastAMG< M, X, PI, A >::usesDirectCoarseLevelSolver ( ) const

Check whether the coarse solver used is a direct solver.

Returns
True if the coarse level solver is a direct solver.
template<class M , class X , class PI , class A >
Dune::Amg::FastAMG< M, X, PI, A >::~FastAMG ( )

Variable Documentation

template<class M , class X , class PI = SequentialInformation, class A = std::allocator<X>>
OperatorHierarchy::AggregatesMapList::const_iterator Dune::Amg::FastAMG< M, X, PI, A >::LevelContext::aggregates

The iterator over the aggregates maps.

template<class M , class X , class PI = SequentialInformation, class A = std::allocator<X>>
std::size_t Dune::Amg::FastAMG< M, X, PI, A >::LevelContext::level

The level index.

template<class M , class X , class PI = SequentialInformation, class A = std::allocator<X>>
Hierarchy<Domain,A>::Iterator Dune::Amg::FastAMG< M, X, PI, A >::LevelContext::lhs

The iterator over the left hand side.

template<class M , class X , class PI = SequentialInformation, class A = std::allocator<X>>
OperatorHierarchy::ParallelMatrixHierarchy::ConstIterator Dune::Amg::FastAMG< M, X, PI, A >::LevelContext::matrix

The iterator over the matrices.

template<class M , class X , class PI = SequentialInformation, class A = std::allocator<X>>
ParallelInformationHierarchy::Iterator Dune::Amg::FastAMG< M, X, PI, A >::LevelContext::pinfo

The iterator over the parallel information.

template<class M , class X , class PI = SequentialInformation, class A = std::allocator<X>>
OperatorHierarchy::RedistributeInfoList::const_iterator Dune::Amg::FastAMG< M, X, PI, A >::LevelContext::redist

The iterator over the redistribution information.

template<class M , class X , class PI = SequentialInformation, class A = std::allocator<X>>
Hierarchy<Domain,A>::Iterator Dune::Amg::FastAMG< M, X, PI, A >::LevelContext::residual

The iterator over the residuals.

template<class M , class X , class PI = SequentialInformation, class A = std::allocator<X>>
Hierarchy<Range,A>::Iterator Dune::Amg::FastAMG< M, X, PI, A >::LevelContext::rhs

The iterator over the right hand sided.