petsc-3.3-p7 2013-05-11

PCLSC

Preconditioning for Schur complements, based on Least Squares Commutators

Options Database Key

-pc_lsc_scale_diag -Use the diagonal of A for scaling

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 inv(A00) A01

PCLSC currently doesn't do anything with A11, so let's assume it is 0. The idea is that a good approximation to inv(S) is given by

   inv(A10 A01) A10 A00 A01 inv(A10 A01)

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 A10 is a Laplacian, call it L. The current interface is to hang L and a preconditioning matrix Lp on the preconditioning matrix.

If you had called KSPSetOperators(ksp,S,Sp,flg), 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

   PetscObjectCompose((PetscObject)Sp,"LSC_L",(PetscObject)L);
   PetscObjectCompose((PetscObject)Sp,"LSC_Lp",(PetscObject)Lp);

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. ML's 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

See Also

PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCFIELDSPLIT,
PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition(), MatCreateSchurComplement()

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

Examples

src/snes/examples/tutorials/ex55.c.html
Elman, Howle, Shadid, Shuttleworth, and Tuminaro, Block preconditioners based on approximate commutators, 2006.- - Silvester, Elman, Kay, Wathen, Efficient preconditioning of the linearized Navier-Stokes equations for incompressible flow, 2001.