PCFieldSplitSetSchurFactType#

sets which blocks of the approximate block factorization to retain in the preconditioner [MGW00] and [Ips01]

Synopsis#

#include "petscpc.h" 
PetscErrorCode PCFieldSplitSetSchurFactType(PC pc, PCFieldSplitSchurFactType ftype)

Collective

Input Parameters#

Options Database Key#

  • -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> - default is full

Notes#

The FULL factorization is

(ABCE)=(10CA1I)(A00S)(IA1B0I)=LDU.\left(\begin{array}{cc} A & B \\ C & E \\ \end{array}\right) = \left(\begin{array}{cc} 1 & 0 \\ C*A^{-1} & I \\ \end{array}\right) \left(\begin{array}{cc} A & 0 \\ 0 & S \\ \end{array}\right) \left(\begin{array}{cc} I & A^{-1}B \\ 0 & I \\ \end{array}\right) = L D U.

where S=ECA1B S = E - C*A^{-1}*B . In practice, the full factorization is applied via block triangular solves with the grouping L(DU)L*(D*U). UPPER uses DUD*U, LOWER uses LDL*D, and DIAG is the diagonal part with the sign of S S flipped (because this makes the preconditioner positive definite for many formulations, thus allowing the use of KSPMINRES). Sign flipping of S S can be turned off with PCFieldSplitSetSchurScale().

If AA and SS are solved exactly

  • 1 - FULL factorization is a direct solver.

  • 2 - The preconditioned operator with LOWER or UPPER has all eigenvalues equal to 1 and minimal polynomial of degree 2, so KSPGMRES converges in 2 iterations.

  • 3 - With DIAG, the preconditioned operator has three distinct nonzero eigenvalues and minimal polynomial of degree at most 4, so KSPGMRES converges in at most 4 iterations.

If the iteration count is very low, consider using KSPFGMRES or KSPGCR which can use one less preconditioner application in this case. Note that the preconditioned operator may be highly non-normal, so such fast convergence may not be observed in practice.

For symmetric problems in which AA is positive definite and SS is negative definite, DIAG can be used with KSPMINRES.

A flexible method like KSPFGMRES or KSPGCR, Flexible Krylov Methods, must be used if the fieldsplit preconditioner is nonlinear (e.g. a few iterations of a Krylov method is used to solve with A or S).

References#

[Ips01]

Ilse CF Ipsen. A note on preconditioning nonsymmetric matrices. SIAM Journal on Scientific Computing, 23(3):1050–1051, 2001.

[MGW00]

Malcolm F Murphy, Gene H Golub, and Andrew J Wathen. A note on preconditioning for indefinite linear systems. SIAM Journal on Scientific Computing, 21(6):1969–1972, 2000.

See Also#

Solving Block Matrices with PCFIELDSPLIT, PC, PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCFieldSplitSetSchurScale(), Flexible Krylov Methods, PCFieldSplitSetSchurPre()

Level#

intermediate

Location#

src/ksp/pc/impls/fieldsplit/fieldsplit.c

Examples#

src/dm/impls/stag/tutorials/ex4.c

Implementations#

PCFieldSplitSetSchurFactType_FieldSplit() in src/ksp/pc/impls/fieldsplit/fieldsplit.c


Index of all PC routines
Table of Contents for all manual pages
Index of all manual pages