petsc-3.6.1 2015-08-06
Report Typos and Errors

PCFieldSplitSetSchurPre

Indicates if the Schur complement is preconditioned by a preconditioner constructed by the A11 matrix. Otherwise no preconditioner is used.

Synopsis

#include "petscpc.h" 
PetscErrorCode PCFieldSplitSetSchurPre(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
Collective on PC

Input Parameters

pc - the preconditioner context
ptype - which matrix to use for preconditioning the Schur complement: PC_FIELDSPLIT_SCHUR_PRE_A11 (default), PC_FIELDSPLIT_SCHUR_PRE_SELF, PC_FIELDSPLIT_PRE_USER
userpre - matrix to use for preconditioning, or NULL

Options Database

-pc_fieldsplit_schur_precondition <self,selfp,user,a11,full> default is a11 - Notes:
   If ptype is
       user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument
            to this function).
       a11 then the preconditioner for the Schur complement is generated by the block diagonal part of the preconditioner
            matrix associated with the Schur complement (i.e. A11), not he Schur complement matrix
       full then the preconditioner uses the exact Schur complement (this is expensive)
       self the preconditioner for the Schur complement is generated from the Schur complement matrix itself:
            The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC
            preconditioner
       selfp then the preconditioning matrix is an explicitly-assembled approximation Sp = A11 - A10 inv(diag(A00)) A01
            This is only a good preconditioner when diag(A00) is a good preconditioner for A00. Optionally, A00 can be
            lumped before extracting the diagonal: -fieldsplit_1_mat_schur_complement_ainv_type lump; diag is the default.

When solving a saddle point problem, where the A11 block is identically zero, using a11 as the ptype only makes sense with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and -fieldsplit_1_pc_type lsc which uses the least squares commutator to compute a preconditioner for the Schur complement.

See Also

PCFieldSplitGetSchurPre(), PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType,
MatSchurComplementSetAinvType(), PCLSC

Level:intermediate
Location:
src/ksp/pc/impls/fieldsplit/fieldsplit.c
Index of all PC routines
Table of Contents for all manual pages
Index of all manual pages

Examples

src/snes/examples/tutorials/ex70.c.html