| KASKADE 7 development version
    | 
An inexact Uzawa type solver for solving symmetric saddle point systems. More...
#include <uzawa.hh>
An inexact Uzawa type solver for solving symmetric saddle point systems.
\[ \begin{array}{cc} A & B^T \\ B & \end{array}\begin{matrix}{c} u \\ \lambda \end{matrix} = \begin{matrix}{c} f \\ g \end{matrix} \]
It implements an inexact preconditioned conjugate gradient method for the Schur complement \( S = BA^{-1}B^T\). The iteration is terminated once the preconditioned residual energy drops below a prescribed reduction factor times the solution energy, i.e. \( r^T\hat S^{-1} S \hat S^{-1} r \le \mathrm{TOL}^2 p^T S p \).
| X | the linear space of the primal variable \( u \) | 
| Y | the linear space of the dual variable \( \lambda \) | 
 
 | Public Types | |
| using | domain_type = LinearProductSpace< typename X::field_type, boost::fusion::vector< X, Y > > | 
| The domain type \( X \times Y \) of the Uzawa iteration.  More... | |
| using | range_type = domain_type | 
| typedef X::field_type | field_type | 
| typedef field_type | Scalar | 
| typedef domain_type | Domain | 
| typedef range_type | Range | 
| Public Member Functions | |
| UzawaSolver (Dune::LinearOperator< X, X > &A_, Dune::InverseOperator< X, X > &invA_, Dune::LinearOperator< X, Y > &B_, Dune::LinearOperator< Y, X > &Bt_, Dune::Preconditioner< Y, Y > &precS_, double reduction_, int maxit_, int verbose_) | |
| Constructor.  More... | |
| virtual void | apply (Domain &x, Range &b, Dune::InverseOperatorResult &res) override | 
| Solves the saddle point system \( K x = b \).  More... | |
| virtual void | apply (Domain &x, Range &b, double reduction_, Dune::InverseOperatorResult &res) | 
| void | apply (Domain &x, Range &b) | 
| virtual Dune::SolverCategory::Category | category () const override | 
| returns the category of the operator  More... | |
| typedef domain_type Kaskade::UzawaSolver< X, Y >::Domain | 
| using Kaskade::UzawaSolver< X, Y >::domain_type = LinearProductSpace<typename X::field_type,boost::fusion::vector<X,Y> > | 
The domain type \( X \times Y \) of the Uzawa iteration.
This is the LinearProductSpace of X and Y
| typedef X::field_type Kaskade::UzawaSolver< X, Y >::field_type | 
| typedef range_type Kaskade::UzawaSolver< X, Y >::Range | 
| using Kaskade::UzawaSolver< X, Y >::range_type = domain_type | 
| typedef field_type Kaskade::UzawaSolver< X, Y >::Scalar | 
| 
 | inline | 
Constructor.
| A_ | the linear operator \( A : X \to X \) | 
| invA_ | the (approximate) inverse operator \( A^{-1}: X \to X \), needs to be a subclass of Dune::InverseOperator<X,X> | 
| B_ | the linear operaotr \( B: X\to Y \) | 
| Bt_ | the linear operator \( B^T : Y\to X\) | 
| precS_ | the Schur complement preconditioner \( P: Y\to Y \), \( P\approx S \) | 
| reduction_ | the default required reduction in preconditioned Schur complement's residual | 
| maxit_ | the maximum number of iterations | 
| verbose_ | printing level to stdout (0=silent, 1=result, 2=every iteration) | 
All the given linear operators and the Schur preconditioner are kept by reference, and have to exist during the lifetime of the UzawaSolver object.
| 
 | inline | 
| 
 | inlinevirtual | 
| 
 | inlineoverridevirtual | 
Solves the saddle point system \( K x = b \).
| [in,out] | x | solution. On entry, an initial guess is provided. | 
| [in] | b | right hand side. This will not be modified. | 
Definition at line 112 of file uzawa.hh.
Referenced by Kaskade::UzawaSolver< X, Y >::apply().
| 
 | inlineoverridevirtual |