:orphan:
# MatILUFactor
Performs in-place ILU factorization of matrix. 
## Synopsis
```
#include "petscmat.h" 
PetscErrorCode MatILUFactor(Mat mat, IS row, IS col, const MatFactorInfo *info)
```
Collective


## Input Parameters

- ***mat -*** the matrix
- ***row -*** row permutation
- ***col -*** column permutation
- ***info -*** structure containing


```none
      levels - number of levels of fill.
      expected fill - as ratio of original fill.
      1 or 0 - indicating force fill on diagonal (improves robustness for matrices
                missing diagonal entries)
```





## Notes
Most users should employ the `KSP` interface for linear solvers
instead of working directly with matrix algebra routines such as this.
See, e.g., `KSPCreate()`.

Probably really in-place only when level of fill is zero, otherwise allocates
new space to store factored matrix and deletes previous memory. The preferred approach is to use `MatGetFactor()`, `MatILUFactorSymbolic()`, and `MatILUFactorNumeric()`
when not using `KSP`.


## Developer Note
The Fortran interface is not autogenerated as the
interface definition cannot be generated correctly [due to MatFactorInfo]


## See Also
 [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `MatILUFactorSymbolic()`, `MatLUFactorNumeric()`, `MatCholeskyFactor()`, `MatFactorInfo`

## Level
developer

## Location
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/mat/interface/matrix.c.html#MatILUFactor">src/mat/interface/matrix.c</A>

## Implementations

<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/mat/impls/aij/seq/aij.c.html#MatILUFactor_SeqAIJ">MatILUFactor_SeqAIJ in src/mat/impls/aij/seq/aij.c</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/mat/impls/baij/seq/baij.c.html#MatILUFactor_SeqBAIJ">MatILUFactor_SeqBAIJ in src/mat/impls/baij/seq/baij.c</A><BR>


---
[Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/mat/interface/matrix.c)


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