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#
pc - the preconditioner context
ftype - which blocks of factorization to retain,
PC_FIELDSPLIT_SCHUR_FACT_FULL
is default
Options Database Key#
-pc_fieldsplit_schur_fact_type <diag,lower,upper,full> - default is
full
Notes#
The FULL factorization is
(A B) = (1 0) (A 0) (1 Ainv*B) = L D U
(C E) (C*Ainv 1) (0 S) (0 1)
.vb
where S = E - C*Ainv*B. In practice, the full factorization is applied via block triangular solves with the grouping $L*(D*U)$. UPPER uses $D*U$, LOWER uses $L*D$,
and DIAG is the diagonal part with the sign of S flipped (because this makes the preconditioner positive definite for many formulations,
thus allowing the use of `KSPMINRES)`. Sign flipping of S can be turned off with `PCFieldSplitSetSchurScale()`.
If A and S are solved exactly
.vb
*) FULL factorization is a direct solver.
*) 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.
*) 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 A is positive definite and S 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, PC
, PCFieldSplitGetSubKSP()
, PCFIELDSPLIT
, PCFieldSplitSetFields()
, PCFieldSplitSchurPreType
, PCFieldSplitSetSchurScale()
,
Flexible Krylov Methods
Level#
intermediate
Location#
Examples#
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