:orphan:
# PCFieldSplitSetSchurPre
Indicates from what operator the preconditioner is constructucted for the Schur complement. The default is the A11 matrix. 
## Synopsis
```
#include "petscpc.h" 
PetscErrorCode PCFieldSplitSetSchurPre(PC pc, PCFieldSplitSchurPreType ptype, Mat pre)
```
Collective


## Input Parameters

- ***pc      -*** the preconditioner context
- ***ptype   -*** which matrix to use for preconditioning the Schur complement: `PC_FIELDSPLIT_SCHUR_PRE_A11` (default),
`PC_FIELDSPLIT_SCHUR_PRE_SELF`, `PC_FIELDSPLIT_SCHUR_PRE_USER`,
`PC_FIELDSPLIT_SCHUR_PRE_SELFP`, and `PC_FIELDSPLIT_SCHUR_PRE_FULL`
- ***pre -*** matrix to use for preconditioning, or `NULL`



## Options Database Keys

- ***-pc_fieldsplit_schur_precondition <self,selfp,user,a11,full> -*** default is a11. See notes for meaning of various arguments
- ***-fieldsplit_1_pc_type <pctype> -*** the preconditioner algorithm that is used to construct the preconditioner from the operator





## Notes
If ptype is

- ***a11 -*** the preconditioner for the Schur complement is generated from the block diagonal part of the preconditioner
matrix associated with the Schur complement (i.e. A11), not the Schur complement matrix
- ***self -*** the preconditioner for the Schur complement is generated from the symbolic representation of the Schur complement matrix:
The only preconditioner that currently works with this symbolic representation matrix object is the `PCLSC`
preconditioner
- ***user -*** the preconditioner for the Schur complement is generated from the user provided matrix (pre argument
to this function).
- ***selfp -*** the preconditioning for the Schur complement is generated from an explicitly-assembled approximation Sp = A11 - A10 inv(diag(A00)) A01
This is only a good preconditioner when diag(A00) is a good preconditioner for A00. Optionally, A00 can be
lumped before extracting the diagonal using the additional option -fieldsplit_1_mat_schur_complement_ainv_type lump
- ***full -*** the preconditioner for the Schur complement is generated from the exact Schur complement matrix representation
computed internally by `PCFIELDSPLIT` (this is expensive)
useful mostly as a test that the Schur complement approach can work for your problem


When solving a saddle point problem, where the A11 block is identically zero, using a11 as the ptype only makes sense
with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
-fieldsplit_1_pc_type lsc which uses the least squares commutator to compute a preconditioner for the Schur complement.


## See Also
 [](sec_block_matrices), `PC`, `PCFieldSplitGetSchurPre()`, `PCFieldSplitGetSubKSP()`, `PCFIELDSPLIT`, `PCFieldSplitSetFields()`, `PCFieldSplitSchurPreType`,
`MatSchurComplementSetAinvType()`, `PCLSC`,
`PCFieldSplitSchurPreType`

## Level
intermediate

## Location
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/pc/impls/fieldsplit/fieldsplit.c.html#PCFieldSplitSetSchurPre">src/ksp/pc/impls/fieldsplit/fieldsplit.c</A>

## Examples
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/dm/impls/stag/tutorials/ex4.c.html">src/dm/impls/stag/tutorials/ex4.c</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/snes/tutorials/ex70.c.html">src/snes/tutorials/ex70.c</A><BR>

## Implementations

<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/pc/impls/fieldsplit/fieldsplit.c.html#PCFieldSplitSetSchurPre_FieldSplit">PCFieldSplitSetSchurPre_FieldSplit in src/ksp/pc/impls/fieldsplit/fieldsplit.c</A><BR>


---
[Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/ksp/pc/impls/fieldsplit/fieldsplit.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)  
