SYNOPSIS

    template <class Matrix, class Vector, class Solver, class Preconditioner, class Size, class Real>
    int pcg_abtb (const Matrix& A, const Matrix& B, Vector& u, Vector& p,
      const Vector& Mf, const Vector& Mg, const Preconditioner& S1,
      const Solver& inner_solver_A, Size& max_iter, Real& tol,
      odiststream *p_derr = 0, std::string label = "pcg_abtb");

    template <class Matrix, class Vector, class Solver, class Preconditioner, class Size, class Real>
    int pcg_abtbc (const Matrix& A, const Matrix& B, const Matrix& C, Vector& u, Vector& p,
      const Vector& Mf, const Vector& Mg, const Preconditioner& S1,
      const Solver& inner_solver_A, Size& max_iter, Real& tol,
      odiststream *p_derr = 0, std::string label = "pcg_abtbc");

The synopsis is the same with the pminres algorithm.

EXAMPLES

See the user's manual for practical examples for the nearly incompressible elasticity, the Stokes and the Navier-Stokes problems.

DESCRIPTION

Preconditioned conjugate gradient algorithm on the pressure p applied to the stabilized stokes problem:

       [ A  B^T ] [ u ]    [ Mf ]
       [        ] [   ]  = [    ]
       [ B  -C  ] [ p ]    [ Mg ]

where A is symmetric positive definite and C is symmetric positive and semi-definite. Such mixed linear problems appears for instance with the discretization of Stokes problems with stabilized P1-P1 element, or with nearly incompressible elasticity. Formaly u = inv(A)*(Mf - B^T*p) and the reduced system writes for all non-singular matrix S1:

     inv(S1)*(B*inv(A)*B^T)*p = inv(S1)*(B*inv(A)*Mf - Mg)

Uzawa or conjugate gradient algorithms are considered on the reduced problem. Here, S1 is some preconditioner for the Schur complement S=B*inv(A)*B^T. Both direct or iterative solvers for S1*q = t are supported. Application of inv(A) is performed via a call to a solver for systems such as A*v = b. This last system may be solved either by direct or iterative algorithms, thus, a general matrix solver class is submitted to the algorithm. For most applications, such as the Stokes problem, the mass matrix for the p variable is a good S1 preconditioner for the Schur complement. The stoping criteria is expressed using the S1 matrix, i.e. in L2 norm when this choice is considered. It is scaled by the L2 norm of the right-hand side of the reduced system, also in S1 norm.