# MatGetSchurComplement
Obtain the Schur complement from eliminating part of the matrix in another part. 
## Synopsis
```
#include "petscksp.h" 
PetscErrorCode MatGetSchurComplement(Mat A, IS isrow0, IS iscol0, IS isrow1, IS iscol1, MatReuse mreuse, Mat *S, MatSchurComplementAinvType ainvtype, MatReuse preuse, Mat *Sp)
```
Collective


## Input Parameters

- ***A      -*** matrix in which the complement is to be taken
- ***isrow0 -*** rows to eliminate
- ***iscol0 -*** columns to eliminate, (isrow0,iscol0) should be square and nonsingular
- ***isrow1 -*** rows in which the Schur complement is formed
- ***iscol1 -*** columns in which the Schur complement is formed
- ***mreuse -*** `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`, use `MAT_IGNORE_MATRIX` to put nothing in S
- ***ainvtype -*** the type of approximation used for the inverse of the (0,0) block used in forming Sp:
`MAT_SCHUR_COMPLEMENT_AINV_DIAG`, `MAT_SCHUR_COMPLEMENT_AINV_LUMP`, `MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG`, or `MAT_SCHUR_COMPLEMENT_AINV_FULL`
- ***preuse -*** `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`, use `MAT_IGNORE_MATRIX` to put nothing in Sp



## Output Parameters

- ***S      -*** exact Schur complement, often of type `MATSCHURCOMPLEMENT` which is difficult to use for preconditioning
- ***Sp     -*** approximate Schur complement from which a preconditioner can be built A11 - A10 inv(DIAGFORM(A00)) A01





## Notes
Since the real Schur complement is usually dense, providing a good approximation to Sp usually requires
application-specific information.

Sometimes users would like to provide problem-specific data in the Schur complement, usually only for special row
and column index sets.  In that case, the user should call `PetscObjectComposeFunction()` on the *S matrix and pass mreuse of `MAT_REUSE_MATRIX` to set
"MatGetSchurComplement_C" to their function.  If their function needs to fall back to the default implementation, it
should call `MatGetSchurComplement_Basic()`.

`MatCreateSchurComplement()` takes as arguments the four submatrices and returns the virtual Schur complement (what this function returns in S).

`MatSchurComplementGetPmat()` takes the virtual Schur complement and returns an explicit approximate Schur complement (what this returns in Sp).

In other words calling `MatCreateSchurComplement()` followed by `MatSchurComplementGetPmat()` produces the same output as this function but with slightly different
inputs. The actually submatrices of the original block matrix instead of index sets to the submatrices.


## Developer Note
The API that includes `MatGetSchurComplement()`, `MatCreateSchurComplement()`, `MatSchurComplementGetPmat()` should be refactored to
remove redundancy and be clearer and simpler.


## See Also
 [](chapter_ksp), `MatCreateSubMatrix()`, `PCFIELDSPLIT`, `MatCreateSchurComplement()`, `MatSchurComplementAinvType`

## Level
advanced

## Location
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/utils/schurm/schurm.c.html#MatGetSchurComplement">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](/docs/manualpages/index.md)  
[Index of all manual pages](/docs/manualpages/singleindex.md)  
