| FEAT 3
    Finite Element Analysis Toolbox | 
Uzawa preconditioner. More...
#include <uzawa_precond.hpp>
 
  
| Public Types | |
| typedef SolverBase< VectorType > | BaseClass | 
| base-class typedef  More... | |
| typedef MatrixA_::DataType | DataType | 
| our data type  More... | |
| typedef SolverBase< VectorTypeV > | SolverA | 
| A-block solver type.  More... | |
| typedef SolverBase< VectorTypeP > | SolverS | 
| S-block solver type.  More... | |
| typedef LAFEM::TupleVector< VectorTypeV, VectorTypeP > | VectorType | 
| system vector type  More... | |
| typedef MatrixD_::VectorTypeL | VectorTypeP | 
| pressure vector type  More... | |
| typedef MatrixB_::VectorTypeL | VectorTypeV | 
| velocity vector type  More... | |
| Public Member Functions | |
| UzawaPrecond (const MatrixA_ &matrix_a, const MatrixB_ &matrix_b, const MatrixD_ &matrix_d, const FilterV_ &filter_v, const FilterP_ &filter_p, std::shared_ptr< SolverA > solver_a, std::shared_ptr< SolverS > solver_s, UzawaType type=UzawaType::diagonal, bool auto_init_s=true) | |
| Constructs a Uzawa preconditioner.  More... | |
| virtual Status | apply (LAFEM::TupleVector< MatrixB_::VectorTypeL, MatrixD_::VectorTypeL > &vec_cor, const LAFEM::TupleVector< MatrixB_::VectorTypeL, MatrixD_::VectorTypeL > &vec_def)=0 | 
| Solver application method.  More... | |
| virtual Status | apply (VectorType &vec_cor, const VectorType &vec_def) override | 
| virtual void | done () | 
| Finalization method.  More... | |
| virtual void | done_numeric () override | 
| Numeric finalization method.  More... | |
| virtual void | done_symbolic () override | 
| Symbolic finalization method.  More... | |
| virtual void | init () | 
| Initialization method.  More... | |
| virtual void | init_numeric () override | 
| Numeric initialization method.  More... | |
| virtual void | init_symbolic () override | 
| Symbolic initialization method.  More... | |
| virtual String | name () const override | 
| Returns a descriptive string.  More... | |
| Private Attributes | |
| const bool | _auto_init_s | 
| auto-initialize of S-solver  More... | |
| const FilterP_ & | _filter_p | 
| const FilterV_ & | _filter_v | 
| our system filter  More... | |
| const MatrixA_ & | _matrix_a | 
| const MatrixB_ & | _matrix_b | 
| const MatrixD_ & | _matrix_d | 
| std::shared_ptr< SolverA > | _solver_a | 
| our A-block solver  More... | |
| std::shared_ptr< SolverS > | _solver_s | 
| our S-block solver  More... | |
| UzawaType | _uzawa_type | 
| our Uzawa type  More... | |
| VectorTypeP | _vec_tmp_p | 
| VectorTypeV | _vec_tmp_v | 
| a temporary defect vector  More... | |
Uzawa preconditioner.
This class implements the Uzawa preconditioner, which is a special preconditioner for saddle-point systems of the form
\[\begin{bmatrix} A & B\\D & 0\end{bmatrix} \cdot \begin{bmatrix} x_u\\x_p\end{bmatrix} = \begin{bmatrix} f_u\\f_p\end{bmatrix}\]
Let \( S \approx -DA^{-1}B\) be an approximation of the Schur-complement matrix, then this class supports a total of four different preconditioners for the saddle-point system above:
\[\begin{bmatrix} A & 0\\0 & S\end{bmatrix} \cdot \begin{bmatrix} x_u\\x_p\end{bmatrix} = \begin{bmatrix} f_u\\f_p\end{bmatrix}\]
\[\begin{bmatrix} A & 0\\D & S\end{bmatrix} \cdot \begin{bmatrix} x_u\\x_p\end{bmatrix} = \begin{bmatrix} f_u\\f_p\end{bmatrix}\]
\[\begin{bmatrix} A & B\\0 & S\end{bmatrix} \cdot \begin{bmatrix} x_u\\x_p\end{bmatrix} = \begin{bmatrix} f_u\\f_p\end{bmatrix}\]
\[\begin{bmatrix} I & 0\\-DA^{-1} & I\end{bmatrix} \cdot \begin{bmatrix} A & B\\0 & S\end{bmatrix} \cdot \begin{bmatrix} x_u\\x_p\end{bmatrix} = \begin{bmatrix} f_u\\f_p\end{bmatrix}\]
The required solution steps of \( A^{-1} \) and \( S^{-1} \) are performed by two sub-solvers, which have to be created by the user and supplied to the constructor of this object.
| MatrixA_,MatrixB_,MatrixD_ | The types of the sub-matrices A, B and D. | 
| FilterV_,FilterP_ | The types of the filters for the velocity and pressure components, respectively. | 
Definition at line 76 of file uzawa_precond.hpp.
| typedef SolverBase<VectorType> FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::BaseClass | 
base-class typedef
Definition at line 89 of file uzawa_precond.hpp.
| typedef MatrixA_::DataType FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::DataType | 
our data type
Definition at line 81 of file uzawa_precond.hpp.
| typedef SolverBase<VectorTypeV> FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::SolverA | 
A-block solver type.
Definition at line 92 of file uzawa_precond.hpp.
| typedef SolverBase<VectorTypeP> FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::SolverS | 
S-block solver type.
Definition at line 94 of file uzawa_precond.hpp.
| typedef LAFEM::TupleVector<VectorTypeV, VectorTypeP> FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::VectorType | 
system vector type
Definition at line 87 of file uzawa_precond.hpp.
| typedef MatrixD_::VectorTypeL FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::VectorTypeP | 
pressure vector type
Definition at line 85 of file uzawa_precond.hpp.
| typedef MatrixB_::VectorTypeL FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::VectorTypeV | 
velocity vector type
Definition at line 83 of file uzawa_precond.hpp.
| 
 | inlineexplicit | 
Constructs a Uzawa preconditioner.
| [in] | matrix_a,matrix_b,matrix_d | The three sub-matrices of the saddle-point matrix. | 
| [in] | filter_v,filter_p | The velocity and pressure filter, resp. | 
| [in] | solver_a | The solver representing \(A^{-1}\). | 
| [in] | solver_s | The solver representing \(S^{-1}\). | 
| [in] | type | Specifies the type of the preconditioner. See this class' documentation for details. | 
| [in] | auto_init_s | Specifies whether this solver object should call the init/done functions of the solver_sobject. If set tofalse, then the caller is responsible for the initialization and finalization of the S-block solver object. | 
Definition at line 139 of file uzawa_precond.hpp.
References XASSERTM.
| 
 | inlinevirtual | 
Definition at line 163 of file uzawa_precond.hpp.
| 
 | pure virtualinherited | 
Solver application method.
This method applies the solver represented by this object onto a given defect vector and returns the corresponding correction vector.
correct() method which corrects an initial solution instead of starting with the null vector.| [out] | vec_cor | The vector that shall receive the solution of the linear system. It is assumed to be allocated, but its numerical contents may be undefined upon calling this method. | 
| [in] | vec_def | The vector that represents the right-hand-side of the linear system to be solved. | 
| 
 | inlineoverridevirtual | 
Definition at line 224 of file uzawa_precond.hpp.
| 
 | inlinevirtualinherited | 
| 
 | inlineoverridevirtual | 
Numeric finalization method.
This method is called to release any data allocated in the numeric initialization step.
Reimplemented from FEAT::Solver::SolverBase< LAFEM::TupleVector< MatrixB_::VectorTypeL, MatrixD_::VectorTypeL > >.
Definition at line 199 of file uzawa_precond.hpp.
References FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::_auto_init_s, FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::_solver_a, FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::_solver_s, and FEAT::Solver::SolverBase< Vector_ >::done_numeric().
| 
 | inlineoverridevirtual | 
Symbolic finalization method.
This method is called to release any data allocated in the symbolic initialization step.
Reimplemented from FEAT::Solver::SolverBase< LAFEM::TupleVector< MatrixB_::VectorTypeL, MatrixD_::VectorTypeL > >.
Definition at line 209 of file uzawa_precond.hpp.
References FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::_auto_init_s, FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::_solver_a, FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::_solver_s, FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::_uzawa_type, FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::_vec_tmp_v, FEAT::Solver::diagonal, and FEAT::Solver::SolverBase< Vector_ >::done_symbolic().
| 
 | inlinevirtualinherited | 
| 
 | inlineoverridevirtual | 
Numeric initialization method.
This method is called to perform numeric initialization of the solver.
Before this function can be called, the symbolic initialization must be performed. 
Reimplemented from FEAT::Solver::SolverBase< LAFEM::TupleVector< MatrixB_::VectorTypeL, MatrixD_::VectorTypeL > >.
Definition at line 189 of file uzawa_precond.hpp.
References FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::_auto_init_s, FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::_solver_a, FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::_solver_s, and FEAT::Solver::SolverBase< Vector_ >::init_numeric().
| 
 | inlineoverridevirtual | 
Symbolic initialization method.
This method is called to perform symbolic initialization of the solver.
Reimplemented from FEAT::Solver::SolverBase< LAFEM::TupleVector< MatrixB_::VectorTypeL, MatrixD_::VectorTypeL > >.
Definition at line 172 of file uzawa_precond.hpp.
References FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::_auto_init_s, FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::_solver_a, FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::_solver_s, FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::_uzawa_type, FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::_vec_tmp_v, FEAT::Solver::diagonal, and FEAT::Solver::SolverBase< Vector_ >::init_symbolic().
| 
 | inlineoverridevirtual | 
Returns a descriptive string.
Implements FEAT::Solver::SolverBase< LAFEM::TupleVector< MatrixB_::VectorTypeL, MatrixD_::VectorTypeL > >.
Definition at line 167 of file uzawa_precond.hpp.
| 
 | private | 
auto-initialize of S-solver
Definition at line 110 of file uzawa_precond.hpp.
Referenced by FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::done_numeric(), FEAT::Solver::UzawaPrecond< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >::done_numeric(), FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::done_symbolic(), FEAT::Solver::UzawaPrecond< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >::done_symbolic(), FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::init_numeric(), FEAT::Solver::UzawaPrecond< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >::init_numeric(), FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::init_symbolic(), and FEAT::Solver::UzawaPrecond< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >::init_symbolic().
| 
 | private | 
Definition at line 102 of file uzawa_precond.hpp.
| 
 | private | 
our system filter
Definition at line 101 of file uzawa_precond.hpp.
| 
 | private | 
Definition at line 97 of file uzawa_precond.hpp.
| 
 | private | 
Definition at line 98 of file uzawa_precond.hpp.
| 
 | private | 
Definition at line 99 of file uzawa_precond.hpp.
| 
 | private | 
our A-block solver
Definition at line 104 of file uzawa_precond.hpp.
Referenced by FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::done_numeric(), FEAT::Solver::UzawaPrecond< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >::done_numeric(), FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::done_symbolic(), FEAT::Solver::UzawaPrecond< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >::done_symbolic(), FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::init_numeric(), FEAT::Solver::UzawaPrecond< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >::init_numeric(), FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::init_symbolic(), and FEAT::Solver::UzawaPrecond< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >::init_symbolic().
| 
 | private | 
our S-block solver
Definition at line 106 of file uzawa_precond.hpp.
Referenced by FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::done_numeric(), FEAT::Solver::UzawaPrecond< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >::done_numeric(), FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::done_symbolic(), FEAT::Solver::UzawaPrecond< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >::done_symbolic(), FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::init_numeric(), FEAT::Solver::UzawaPrecond< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >::init_numeric(), FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::init_symbolic(), and FEAT::Solver::UzawaPrecond< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >::init_symbolic().
| 
 | private | 
our Uzawa type
Definition at line 108 of file uzawa_precond.hpp.
Referenced by FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::done_symbolic(), FEAT::Solver::UzawaPrecond< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >::done_symbolic(), FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::init_symbolic(), and FEAT::Solver::UzawaPrecond< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >::init_symbolic().
| 
 | private | 
Definition at line 113 of file uzawa_precond.hpp.
| 
 | private | 
a temporary defect vector
Definition at line 112 of file uzawa_precond.hpp.
Referenced by FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::done_symbolic(), and FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::init_symbolic().