:orphan: # 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 ```none 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 ```none 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 A01 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), 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 ```none PetscObjectCompose((PetscObject)Sp,"LSC_L",(PetscObject)L); PetscObjectCompose((PetscObject)Sp,"LSC_Lp",(PetscObject)Lp); ``` And then your Jacobian assembly would look like ```none 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 ```none -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 - **** -*** 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. ## See Also `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 --- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/ksp/pc/impls/lsc/lsc.c) [Index of all PC routines](index.md) [Table of Contents for all manual pages](/manualpages/index.md) [Index of all manual pages](/manualpages/singleindex.md)