version 3.8
Dumux::SeqUzawa< M, X, Y, l > Class Template Reference

A preconditioner based on the Uzawa algorithm for saddle-point problems of the form \( \begin{pmatrix} A & B \\ C & D \end{pmatrix} \begin{pmatrix} u\\ p \end{pmatrix} = \begin{pmatrix} f\\ g \end{pmatrix} \). More...

#include <dumux/linear/preconditioners.hh>

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

Description

template<class M, class X, class Y, int l = 1>
class Dumux::SeqUzawa< M, X, Y, l >

This preconditioner is especially suited for solving the incompressible (Navier-)Stokes equations. Here, \(D = 0\) and \(B = C^T\) if \(\rho = 1\). We do not expect good convergence if energy or mass transport is considered.

See: Benzi, M., Golub, G. H., & Liesen, J. (2005). Numerical solution of saddle point problems. Acta numerica, 14, 1-137 [14] and
Ho, N., Olson, S. D., & Walker, H. F. (2017). Accelerating the Uzawa algorithm. SIAM Journal on Scientific Computing, 39(5), S461-S476 [38]

Template Parameters
MType of the matrix.
XType of the update.
YType of the defect.
lPreconditioner block level (for compatibility reasons, unused).

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...
 

Public Member Functions

 SeqUzawa (const std::shared_ptr< const Dune::AssembledLinearOperator< M, X, Y > > &op, const Dune::ParameterTree &params)
 Constructor. More...
 
virtual void pre (X &x, Y &b)
 Prepare the preconditioner. More...
 
virtual void apply (X &update, const Y &currentDefect)
 Apply the preconditioner. More...
 
virtual void post (X &x)
 Clean up. More...
 
virtual Dune::SolverCategory::Category category () const
 Category of the preconditioner (see SolverCategory::Category) More...
 

Member Typedef Documentation

◆ domain_type

template<class M , class X , class Y , int l = 1>
using Dumux::SeqUzawa< M, X, Y, l >::domain_type = X

◆ field_type

template<class M , class X , class Y , int l = 1>
using Dumux::SeqUzawa< M, X, Y, l >::field_type = typename X::field_type

◆ matrix_type

template<class M , class X , class Y , int l = 1>
using Dumux::SeqUzawa< M, X, Y, l >::matrix_type = M

◆ range_type

template<class M , class X , class Y , int l = 1>
using Dumux::SeqUzawa< M, X, Y, l >::range_type = Y

◆ scalar_field_type

template<class M , class X , class Y , int l = 1>
using Dumux::SeqUzawa< M, X, Y, l >::scalar_field_type = Dune::Simd::Scalar<field_type>

Constructor & Destructor Documentation

◆ SeqUzawa()

template<class M , class X , class Y , int l = 1>
Dumux::SeqUzawa< M, X, Y, l >::SeqUzawa ( const std::shared_ptr< const Dune::AssembledLinearOperator< M, X, Y > > &  op,
const Dune::ParameterTree &  params 
)
inline
Parameters
opThe linear operator.
paramsCollection of parameters.

Member Function Documentation

◆ apply()

template<class M , class X , class Y , int l = 1>
virtual void Dumux::SeqUzawa< M, X, Y, l >::apply ( X &  update,
const Y &  currentDefect 
)
inlinevirtual
Parameters
updateThe update to be computed.
currentDefectThe current defect.

◆ category()

template<class M , class X , class Y , int l = 1>
virtual Dune::SolverCategory::Category Dumux::SeqUzawa< M, X, Y, l >::category ( ) const
inlinevirtual

◆ post()

template<class M , class X , class Y , int l = 1>
virtual void Dumux::SeqUzawa< M, X, Y, l >::post ( X &  x)
inlinevirtual

◆ pre()

template<class M , class X , class Y , int l = 1>
virtual void Dumux::SeqUzawa< M, X, Y, l >::pre ( X &  x,
Y &  b 
)
inlinevirtual

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