:orphan:
# MatSchurComplementSetSubMatrices
Sets the matrices that define the Schur complement 
## Synopsis
```
#include "petscksp.h" 
PetscErrorCode MatSchurComplementSetSubMatrices(Mat S, Mat A00, Mat Ap00, Mat A01, Mat A10, Mat A11)
```
Collective


## Input Parameters

- ***S                -*** matrix obtained with `MatSetType`(S,`MATSCHURCOMPLEMENT`)
+   A00  - the upper-left block of the original matrix A = [A00 A01; A10 A11]
- ***Ap00 -*** preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A00^{-1}
- ***A01  -*** the upper-right block of the original matrix A = [A00 A01; A10 A11]
- ***A10  -*** the lower-left block of the original matrix A = [A00 A01; A10 A11]
- ***A11  -*** (optional) the lower-right block of the original matrix A = [A00 A01; A10 A11]





## Notes
The Schur complement is NOT explicitly formed! Rather, this
object performs the matrix-vector product of the Schur complement by using formula S = A11 - A10 ksp(A00,Ap00) A01

All four matrices must have the same MPI communicator.

`A00` and `A11` must be square matrices.

This is to be used in the context of code such as
```none
     MatSetType(S,MATSCHURCOMPLEMENT);
     MatSchurComplementSetSubMatrices(S,...);
```

while `MatSchurComplementUpdateSubMatrices()` should only be called after `MatCreateSchurComplement()` or `MatSchurComplementSetSubMatrices()`


## See Also
 [](ch_ksp), `Mat`, `MatCreateNormal()`, `MatMult()`, `MatCreate()`, `MatSchurComplementGetKSP()`, `MatSchurComplementUpdateSubMatrices()`, `MatCreateTranspose()`, `MatCreateSchurComplement()`, `MatGetSchurComplement()`

## Level
intermediate

## Location
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/utils/schurm/schurm.c.html#MatSchurComplementSetSubMatrices">src/ksp/ksp/utils/schurm/schurm.c</A>


---
[Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/ksp/ksp/utils/schurm/schurm.c)


[Index of all KSP routines](index.md)  
[Table of Contents for all manual pages](/manualpages/index.md)  
[Index of all manual pages](/manualpages/singleindex.md)  
