Provided by: librheolef-dev_6.7-6_amd64
NAME
solver_abtb -- direct or iterative solver iterface for mixed linear systems
SYNOPSIS
solver_abtb stokes (a,b,mp); solver_abtb elasticity (a,b,c,mp);
DESCRIPTION
The solver_abtb class provides direct or iterative algorithms for some mixed problem: [ A B^T ] [ u ] [ Mf ] [ ] [ ] = [ ] [ B -C ] [ p ] [ Mg ] where A is symmetric positive definite and C is symmetric positive. By default, iterative algorithms are considered for tridimensional problems and direct methods otherwise. An optional argument can change this behavior. Such mixed linear problems appears for instance with the discretization of Stokes problems. The C matrix can be zero and then the corresponding argument can be omitted when invoking the constructor. Non-zero C matrix appears for of Stokes problems with stabilized P1-P1 element, or for nearly incompressible elasticity problems.
DIRECT ALGORITHM
When the kernel of B^T is not reduced to zero, then the pressure p is defined up to a constant and the system is singular. In the case of iterative methods, this is not a problem. But when using direct method, the system is then completed to impose a constraint on the pressure term and the whole matrix is factored one time for all.
ITERATIVE ALGORITHM
The preconditionned conjugate gradient algorithm is used, where the mp matrix is used as preconditionner. See see mixed_solver(4). The linear sub-systems related to the A matrix are also solved by an iterative algorithm. Use a second optional argument to change this default behavior: a factorization and a direct solver can be considered for these sub- systems.
EXAMPLES
See the user's manual for practical examples for the nearly incompressible elasticity, the Stokes and the Navier-Stokes problems.
IMPLEMENTATION
template <class T, class M = rheo_default_memory_model> class solver_abtb_basic { public: // typedefs: typedef typename csr<T,M>::size_type size_type; // allocators: solver_abtb_basic (); solver_abtb_basic (const csr<T,M>& a, const csr<T,M>& b, const csr<T,M>& mp, const solver_option_type& opt = solver_option_type(), const solver_option_type& sub_opt = solver_option_type()); solver_abtb_basic (const csr<T,M>& a, const csr<T,M>& b, const csr<T,M>& c, const csr<T,M>& mp, const solver_option_type& opt = solver_option_type(), const solver_option_type& sub_opt = solver_option_type()); // accessors: void solve (const vec<T,M>& f, const vec<T,M>& g, vec<T,M>& u, vec<T,M>& p) const; protected: // internal void init(); // data: mutable solver_option_type _opt; mutable solver_option_type _sub_opt; csr<T,M> _a; csr<T,M> _b; csr<T,M> _c; csr<T,M> _mp; solver_basic<T,M> _sA; solver_basic<T,M> _sa; solver_basic<T,M> _smp; bool _need_constraint; }; typedef solver_abtb_basic<Float,rheo_default_memory_model> solver_abtb;
SEE ALSO
mixed_solver(4)