KASKADE 7 development version
Public Member Functions | Static Public Attributes | List of all members
Kaskade::AgglomerationPreconditioner< Operator > Class Template Referenceabstract

An agglomeration-based preconditioner for elliptic problems discretized with higher-order Lagrange finite elements. More...

#include <agglomerationPreconditioner.hh>

Detailed Description

template<class Operator>
class Kaskade::AgglomerationPreconditioner< Operator >

An agglomeration-based preconditioner for elliptic problems discretized with higher-order Lagrange finite elements.

This implements an additive two-level subspace correction preconditioner. The first (fine) level consists of the individual ansatz functions, i.e. a Jacobi preconditioner. The second (coarse) level consists of functions constant one on a single cell, zero on cells not incident to vertices of that cell, and something else in between. A drawback of the approach is that the projected Galerkin matrix on the coarse space is denser than a usual linear FE stiffness matrix. Thus, the coarse space functions span individual subspaces (which realizes a Jacobi smoother on the coarse space), such that we only need to store the diagonal of the projected Galerkin matrix.

For higher order ansatz spaces, the dimension of the coarse space is quite small.

Template Parameters
Operatorthe Galerkin operator

The Galerkin operator type can satisfy either the Dune::AssembledLinearOperator interface or, by specialization, be an AssembledGalerkinOperator of block size 1 x 1. For the Dune::AssembledLinearOperator, the matrix type has to implement the Dune::BCRSMatrix interface, the domain and range types have to implement the Dune::BlockVector interface (e.g. FunctionSpaceElement).

Definition at line 157 of file agglomerationPreconditioner.hh.

Inheritance diagram for Kaskade::AgglomerationPreconditioner< Operator >:
Kaskade::SymmetricPreconditioner< Operator::domain_type, Operator::range_type >

Public Member Functions

template<class Space >
 AgglomerationPreconditioner (Operator const &opA, Space const &space)
 Constructor. More...
 
virtual void apply (domain_type &x, range_type const &y)
 Computes \( x \leftarrow By \). More...
 
virtual field_type applyDp (domain_type &x, range_type const &y)
 Computes \( x \leftarrow By \) and returns \( \langle By, y \rangle \). More...
 
virtual bool requiresInitializedInput () const
 Returns true if the target vector x has to be initialized to zero before calling apply or applyDp. More...
 
virtual void pre (Operator::domain_type &, Operator::range_type &)
 Preconditioner preparation. More...
 
virtual void post (Operator::domain_type &x)
 Preconditioner cleanup. More...
 
virtual field_type applyDp (Operator::domain_type &x, Operator::range_type const &y)=0
 Computes \( x \leftarrow By \) and returns \( \langle By, y \rangle \). More...
 
virtual Dune::SolverCategory::Category category () const override
 returns the category of the operator More...
 

Static Public Attributes

static int const category = Dune::SolverCategory::sequential
 

Constructor & Destructor Documentation

◆ AgglomerationPreconditioner()

template<class Operator >
template<class Space >
Kaskade::AgglomerationPreconditioner< Operator >::AgglomerationPreconditioner ( Operator const &  opA,
Space const &  space 
)
inline

Constructor.

Template Parameters
SpaceLagrangian finite element space
Parameters
opAthe assembled symmetric Galerkin operator
spacethe Lagrangian FE space with which the operator A is discretized

Definition at line 179 of file agglomerationPreconditioner.hh.

Member Function Documentation

◆ apply()

template<class Operator >
virtual void Kaskade::AgglomerationPreconditioner< Operator >::apply ( domain_type &  x,
range_type const &  y 
)
inlinevirtual

Computes \( x \leftarrow By \).

Definition at line 187 of file agglomerationPreconditioner.hh.

◆ applyDp() [1/2]

template<class Operator >
virtual field_type Kaskade::AgglomerationPreconditioner< Operator >::applyDp ( domain_type &  x,
range_type const &  y 
)
inlinevirtual

Computes \( x \leftarrow By \) and returns \( \langle By, y \rangle \).

Definition at line 196 of file agglomerationPreconditioner.hh.

◆ applyDp() [2/2]

virtual field_type Kaskade::SymmetricPreconditioner< Operator::domain_type , Operator::range_type >::applyDp ( Operator::domain_type &  x,
Operator::range_type const &  y 
)
pure virtualinherited

Computes \( x \leftarrow By \) and returns \( \langle By, y \rangle \).

Implemented in Kaskade::JacobiPreconditionerDetail::JacobiPreconditionerBase< Operator::matrix_type, Operator::domain_type, Operator::range_type >.

◆ category()

virtual Dune::SolverCategory::Category Kaskade::SymmetricPreconditioner< Operator::domain_type , Operator::range_type >::category ( ) const
inlineoverridevirtualinherited

returns the category of the operator

From the Dune doxygen documentation it is unclear what this is supposed to mean. We return a dummy here.

Definition at line 242 of file symmetricOperators.hh.

◆ post()

virtual void Kaskade::SymmetricPreconditioner< Operator::domain_type , Operator::range_type >::post ( Operator::domain_type &  x)
inlinevirtualinherited

Preconditioner cleanup.

The provided default implementation does nothing.

Definition at line 225 of file symmetricOperators.hh.

◆ pre()

virtual void Kaskade::SymmetricPreconditioner< Operator::domain_type , Operator::range_type >::pre ( Operator::domain_type &  ,
Operator::range_type &   
)
inlinevirtualinherited

Preconditioner preparation.

The provided default implementation does nothing.

Definition at line 218 of file symmetricOperators.hh.

◆ requiresInitializedInput()

template<class Operator >
virtual bool Kaskade::AgglomerationPreconditioner< Operator >::requiresInitializedInput ( ) const
inlinevirtual

Returns true if the target vector x has to be initialized to zero before calling apply or applyDp.

Implements Kaskade::SymmetricPreconditioner< Operator::domain_type, Operator::range_type >.

Definition at line 205 of file agglomerationPreconditioner.hh.

Member Data Documentation

◆ category

template<class Operator >
int const Kaskade::AgglomerationPreconditioner< Operator >::category = Dune::SolverCategory::sequential
static

Definition at line 170 of file agglomerationPreconditioner.hh.


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