version 3.10-dev
Dumux::Detail::StokesPreconditioner< M, X, Y, l > Class Template Reference

A Stokes preconditioner (saddle-point problem) for the problem \( \begin{pmatrix} A & B \\ C & 0 \end{pmatrix} \begin{pmatrix} u \\ p \end{pmatrix} = \begin{pmatrix} f \\ g \end{pmatrix}, \). More...

#include <dumux/linear/stokes_solver.hh>

Inheritance diagram for Dumux::Detail::StokesPreconditioner< M, X, Y, l >:

Description

template<class M, class X, class Y, int l = 2>
class Dumux::Detail::StokesPreconditioner< M, X, Y, l >

where A is the discrete velocity operator and B and C are discrete gradient and divergence operators. This preconditioner is especially suited for solving the incompressible Stokes equations.

Template Parameters
MType of the matrix.
XType of the update.
YType of the defect.
lPreconditioner block level

Public Types

using matrix_type = M
 The matrix type the preconditioner is for. More...
 
using domain_type = X
 The domain type of the preconditioner. More...
 
using range_type = Y
 The range type of the preconditioner. More...
 
using field_type = typename X::field_type
 The field type of the preconditioner. More...
 
using scalar_field_type = Dune::Simd::Scalar< field_type >
 Scalar type underlying the field_type. More...
 
using PressureLinearOperator = Dune::MatrixAdapter< P, V, V >
 the type of the pressure operator More...
 

Public Member Functions

 StokesPreconditioner (const std::shared_ptr< const Dune::AssembledLinearOperator< M, X, Y > > &fullLinearOperator, const std::shared_ptr< const PressureLinearOperator > &pressureLinearOperator, const Dune::ParameterTree &params)
 Constructor. More...
 
void pre (X &update, Y &currentDefect) override
 Prepare the preconditioner. More...
 
void apply (X &update, const Y &currentDefect) override
 Apply the preconditioner. More...
 
void post (X &update) override
 Clean up. More...
 
Dune::SolverCategory::Category category () const override
 Category of the preconditioner (see SolverCategory::Category) More...
 

Member Typedef Documentation

◆ domain_type

template<class M , class X , class Y , int l = 2>
using Dumux::Detail::StokesPreconditioner< M, X, Y, l >::domain_type = X

◆ field_type

template<class M , class X , class Y , int l = 2>
using Dumux::Detail::StokesPreconditioner< M, X, Y, l >::field_type = typename X::field_type

◆ matrix_type

template<class M , class X , class Y , int l = 2>
using Dumux::Detail::StokesPreconditioner< M, X, Y, l >::matrix_type = M

◆ PressureLinearOperator

template<class M , class X , class Y , int l = 2>
using Dumux::Detail::StokesPreconditioner< M, X, Y, l >::PressureLinearOperator = Dune::MatrixAdapter<P,V,V>

◆ range_type

template<class M , class X , class Y , int l = 2>
using Dumux::Detail::StokesPreconditioner< M, X, Y, l >::range_type = Y

◆ scalar_field_type

template<class M , class X , class Y , int l = 2>
using Dumux::Detail::StokesPreconditioner< M, X, Y, l >::scalar_field_type = Dune::Simd::Scalar<field_type>

Constructor & Destructor Documentation

◆ StokesPreconditioner()

template<class M , class X , class Y , int l = 2>
Dumux::Detail::StokesPreconditioner< M, X, Y, l >::StokesPreconditioner ( const std::shared_ptr< const Dune::AssembledLinearOperator< M, X, Y > > &  fullLinearOperator,
const std::shared_ptr< const PressureLinearOperator > &  pressureLinearOperator,
const Dune::ParameterTree &  params 
)
inline
Parameters
fullLinearOperatorthe Stokes linear operator
pressureLinearOperatorthe linear operator to be used for the pressure space preconditioner
paramsa parameter tree for the preconditioner configuration

Member Function Documentation

◆ apply()

template<class M , class X , class Y , int l = 2>
void Dumux::Detail::StokesPreconditioner< M, X, Y, l >::apply ( X &  update,
const Y &  currentDefect 
)
inlineoverride
Parameters
updateThe update to be computed.
currentDefectThe current defect.

The currentDefect has be be in a consistent representation, Definition 2.3 Blatt and Bastian (2009) https://doi.org/10.1504/IJCSE.2008.021112 The update is initially zero. At exit the update has to be in a consistent representation. This usually requires communication.

◆ category()

template<class M , class X , class Y , int l = 2>
Dune::SolverCategory::Category Dumux::Detail::StokesPreconditioner< M, X, Y, l >::category ( ) const
inlineoverride

◆ post()

template<class M , class X , class Y , int l = 2>
void Dumux::Detail::StokesPreconditioner< M, X, Y, l >::post ( X &  update)
inlineoverride

◆ pre()

template<class M , class X , class Y , int l = 2>
void Dumux::Detail::StokesPreconditioner< M, X, Y, l >::pre ( X &  update,
Y &  currentDefect 
)
inlineoverride

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