Dune Core Modules (unstable)

Parallel algebraic multigrid based on agglomeration. More...

#include <dune/istl/paamg/amg.hh>

Public Types

typedef M Operator
 The matrix operator type.
 
typedef PI 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 > OperatorHierarchy
 The operator hierarchy type.
 
typedef OperatorHierarchy::ParallelInformationHierarchy ParallelInformationHierarchy
 The parallal data distribution hierarchy type.
 
typedef X Domain
 The domain type.
 
typedef X Range
 The range type.
 
typedef InverseOperator< X, X > CoarseSolver
 the type of the coarse solver.
 
typedef S Smoother
 The type of the smoother. More...
 
typedef SmootherTraits< Smoother >::Arguments SmootherArgs
 The argument type for the construction of the smoother.
 
typedef X domain_type
 The domain type of the preconditioner.
 
typedef X range_type
 The range type of the preconditioner.
 
typedef X::field_type field_type
 The field type of the preconditioner.
 

Public Member Functions

 AMG (OperatorHierarchy &matrices, CoarseSolver &coarseSolver, const SmootherArgs &smootherArgs, const Parameters &parms)
 Construct a new amg with a specific coarse solver. More...
 
template<class C >
 AMG (const Operator &fineOperator, const C &criterion, const SmootherArgs &smootherArgs=SmootherArgs(), const ParallelInformation &pinfo=ParallelInformation())
 Construct an AMG with an inexact coarse solver based on the smoother. More...
 
 AMG (std::shared_ptr< const Operator > fineOperator, const ParameterTree &configuration, const ParallelInformation &pinfo=ParallelInformation())
 Constructor an AMG via ParameterTree. More...
 
 AMG (const AMG &amg)
 Copy constructor.
 
void pre (Domain &x, Range &b)
 Prepare the preconditioner. More...
 
void apply (Domain &v, const Range &d)
 Apply one step of the preconditioner to the system A(v)=d. More...
 
virtual SolverCategory::Category category () const
 Category of the preconditioner (see SolverCategory::Category)
 
void post (Domain &x)
 Clean up. More...
 
template<class A1 >
void getCoarsestAggregateNumbers (std::vector< std::size_t, A1 > &cont)
 Get the aggregate number of each unknown on the coarsest level. More...
 
void recalculateHierarchy ()
 Recalculate the matrix hierarchy. More...
 
bool usesDirectCoarseLevelSolver () const
 Check whether the coarse solver used is a direct solver. More...
 
virtual void pre (X &x, X &b)=0
 Prepare the preconditioner. More...
 
virtual void apply (X &v, const X &d)=0
 Apply one step of the preconditioner to the system A(v)=d. More...
 

Detailed Description

template<class M, class X, class S, class PI = SequentialInformation, class A = std::allocator<X>>
class Dune::Amg::AMG< M, X, S, PI, A >

Parallel algebraic multigrid based on agglomeration.

Template Parameters
MThe LinearOperator type which represents the matrix
XThe vector type
SThe smoother type
AAn allocator for X
Todo:
drop the smoother template parameter and replace with dynamic construction

Member Function Documentation

◆ apply()

virtual void Dune::Preconditioner< X, X >::apply ( X &  v,
const X &  d 
)
pure virtualinherited

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.

Implemented in Dune::SeqOverlappingSchwarz< M, X, TM, TD, TA >, Dune::SeqOverlappingSchwarz< M, X, TM, TS, TA >, Dune::SeqOverlappingSchwarz< M, X, TM, TD, TA >, and Dune::SeqOverlappingSchwarz< M, X, TM, TS, TA >.

◆ pre()

virtual void Dune::Preconditioner< X, X >::pre ( X &  x,
X &  b 
)
pure virtualinherited

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.
X x(0.0);
Y b = ...; // rhs
prec.pre(x,b); // prepare the preconditioner
prec.apply(x,b); // can be called multiple times now...
prec.post(x); // cleanup internal state
Base class for matrix free definition of preconditioners.
Definition: preconditioner.hh:33
Parameters
xThe left hand side of the equation.
bThe right hand side of the equation.

Implemented in Dune::SeqOverlappingSchwarz< M, X, TM, TD, TA >, and Dune::SeqOverlappingSchwarz< M, X, TM, TS, TA >.


The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 21, 23:30, 2024)