:orphan:
# MatSchurComplementUpdateSubMatrices
Updates the Schur complement matrix object with new submatrices 
## Synopsis
```
#include "petscksp.h" 
PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat S, Mat A00, Mat Ap00, Mat A01, Mat A10, Mat A11)
```
Collective


## Input Parameters

- ***S                -*** matrix obtained with `MatCreateSchurComplement()` (or `MatSchurSetSubMatrices()`) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
- ***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
All four matrices must have the same MPI communicator

`A00` and  `A11` must be square matrices

All of the matrices provided must have the same sizes as was used with `MatCreateSchurComplement()` or `MatSchurComplementSetSubMatrices()`
though they need not be the same matrices.

This can only be called after `MatCreateSchurComplement()` or `MatSchurComplementSetSubMatrices()`, it cannot be called immediately after `MatSetType`(S,`MATSCHURCOMPLEMENT`);


## Developer Note
This code is almost identical to `MatSchurComplementSetSubMatrices()`. The API should be refactored.


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

## Level
intermediate

## Location
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/utils/schurm/schurm.c.html#MatSchurComplementUpdateSubMatrices">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)  
