PCLSC#

Preconditioning for Schur complements, based on Least Squares Commutators [EHS+06] [SEKW01]

Options Database Key#

  • -pc_lsc_commute - Whether to commute the LSC preconditioner in the style of Olshanskii

  • -pc_lsc_scale_diag - Whether to scale \(BB^T\) products. Will use the inverse of the diagonal of \(Qscale\) or \(A\) if the former is not provided

Notes#

This preconditioner will normally be used with PCFIELDSPLIT to precondition the Schur complement, but it can be used for any Schur complement system. Consider the Schur complement

\[ S = A11 - A10 A00^{-1} A01 \]

PCLSC currently doesn’t do anything with \(A11\), so let’s assume it is 0. The idea is that a good approximation to \(S^{-1}\) is given by

\[ (A10 A01)^{-1} A10 A00 A01 (A10 A01)^{-1} \]

The product \(A10 A01\) can be computed for you, but you can provide it (this is usually more efficient anyway). In the case of incompressible flow, \(A10 A01\) is a Laplacian; call it \(L\). The current interface is to compose \(L\) and a matrix from which to construct its preconditioning \(Lp\) on the preconditioning matrix.

If you had called KSPSetOperators(ksp,S,Sp), \(S\) should have type MATSCHURCOMPLEMENT and \(Sp\) can be any type you like (PCLSC doesn’t use it directly) but should have matrices composed with it, under the names “LSC_L” and “LSC_Lp”. For example, you might have setup code like this

And then your Jacobian assembly would look like

   PetscObjectQuery((PetscObject)Sp,"LSC_L",(PetscObject*)&L);
   PetscObjectQuery((PetscObject)Sp,"LSC_Lp",(PetscObject*)&Lp);
   if (L) { assembly L }
   if (Lp) { assemble Lp }

With this, you should be able to choose LSC preconditioning, using e.g. the PCML algebraic multigrid to solve with L

   -fieldsplit_1_pc_type lsc -fieldsplit_1_lsc_pc_type ml

Since we do not use the values in Sp, you can still put an assembled matrix there to use normal preconditioners.

References#

[EHS+06]

H.C. Elman, V.E. Howle, J. Shadid, R. Shuttleworth, and R. Tuminaro. Block preconditioners based on approximate commutators. SIAM J. Sci. Comput., 27(5):1651–1668, 2006.

[SEKW01]

D. Silvester, H. Elman, D. Kay, and A. Wathen. Efficient preconditioning of the linearized Navier-Stokes equations for incompressible flow. Journal of Computational and Applied Mathematics, 128(1-2):261–279, 2001.

See Also#

KSP: Linear System Solvers, PCCreate(), PCSetType(), PCType, PC, Block_Preconditioners, PCFIELDSPLIT, PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSetSchurPre(), MatCreateSchurComplement(), MatCreateSchurComplement(), MatSchurComplementSetSubMatrices(), MatSchurComplementUpdateSubMatrices(), MatSchurComplementSetAinvType(), MatGetSchurComplement()

Level#

intermediate

Location#

src/ksp/pc/impls/lsc/lsc.c


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