DUNE PDELab (git)

Two grid operator for AMG with Krylov cycle. More...

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

Public Types

typedef AMG::Domain domain_type
 The domain type of the preconditioner.
 
typedef AMG::Range range_type
 The range type of the preconditioner.
 
typedef X::field_type field_type
 The field type of the preconditioner.
 

Public Member Functions

virtual SolverCategory::Category category () const
 Category of the preconditioner (see SolverCategory::Category)
 
 KAmgTwoGrid (AMG &amg, std::shared_ptr< InverseOperator< Domain, Range > > coarseSolver)
 Constructor. More...
 
void pre (typename AMG::Domain &x, typename AMG::Range &b)
 Prepare the preconditioner. More...
 
void post (typename AMG::Domain &x)
 Clean up. More...
 
void apply (typename AMG::Domain &v, const typename AMG::Range &d)
 Apply one step of the preconditioner to the system A(v)=d. More...
 
InverseOperator< Domain, Range > * coarseSolver ()
 Get a pointer to the coarse grid solver. More...
 
void setLevelContext (std::shared_ptr< typename AMG::LevelContext > p)
 Set the level context pointer. More...
 
 ~KAmgTwoGrid ()
 Destructor.
 
virtual void pre (AMG::Domain &x, AMG::Range &b)=0
 Prepare the preconditioner. More...
 
virtual void apply (AMG::Domain &v, const AMG::Range &d)=0
 Apply one step of the preconditioner to the system A(v)=d. More...
 
virtual void post (AMG::Domain &x)=0
 Clean up. More...
 

Detailed Description

template<class AMG>
class Dune::Amg::KAmgTwoGrid< AMG >

Two grid operator for AMG with Krylov cycle.

Template Parameters
AMGThe type of the underlying agglomeration AMG.

Constructor & Destructor Documentation

◆ KAmgTwoGrid()

template<class AMG >
Dune::Amg::KAmgTwoGrid< AMG >::KAmgTwoGrid ( AMG amg,
std::shared_ptr< InverseOperator< Domain, Range > >  coarseSolver 
)
inline

Constructor.

Parameters
amgThe underlying amg. It is used as the storage for the hierarchic data structures.
coarseSolverThe solver used for the coarse grid correction.

Member Function Documentation

◆ apply() [1/2]

virtual void Dune::Preconditioner< AMG::Domain , AMG::Range >::apply ( AMG::Domain &  v,
const AMG::Range &  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.

◆ apply() [2/2]

template<class AMG >
void Dune::Amg::KAmgTwoGrid< AMG >::apply ( typename AMG::Domain v,
const typename AMG::Range d 
)
inline

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.

References Dune::Amg::postsmooth(), and Dune::Amg::presmooth().

◆ coarseSolver()

template<class AMG >
InverseOperator< Domain, Range > * Dune::Amg::KAmgTwoGrid< AMG >::coarseSolver ( )
inline

Get a pointer to the coarse grid solver.

Returns
The coarse grid solver.

◆ post() [1/2]

virtual void Dune::Preconditioner< AMG::Domain , AMG::Range >::post ( AMG::Domain &  x)
pure virtualinherited

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.

◆ post() [2/2]

template<class AMG >
void Dune::Amg::KAmgTwoGrid< AMG >::post ( typename AMG::Domain x)
inline

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.

◆ pre() [1/2]

virtual void Dune::Preconditioner< AMG::Domain , AMG::Range >::pre ( AMG::Domain &  x,
AMG::Range &  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.

◆ pre() [2/2]

template<class AMG >
void Dune::Amg::KAmgTwoGrid< AMG >::pre ( typename AMG::Domain x,
typename AMG::Range b 
)
inline

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
Parameters
xThe left hand side of the equation.
bThe right hand side of the equation.

◆ setLevelContext()

template<class AMG >
void Dune::Amg::KAmgTwoGrid< AMG >::setLevelContext ( std::shared_ptr< typename AMG::LevelContext >  p)
inline

Set the level context pointer.

Parameters
pThe pointer to set it to.

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