:orphan:
# MatLUFactor
Performs in-place LU factorization of matrix.
## Synopsis
```
#include "petscmat.h"
PetscErrorCode MatLUFactor(Mat mat, IS row, IS col, const MatFactorInfo *info)
```
Collective
## Input Parameters
- ***mat -*** the matrix
- ***row -*** row permutation
- ***col -*** column permutation
- ***info -*** options for factorization, includes
```none
fill - expected fill as ratio of original fill.
dtcol - pivot tolerance (0 no pivot, 1 full column pivoting)
Run with the option -info to determine an optimal value to use
```
## 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()`.
This changes the state of the matrix to a factored matrix; it cannot be used
for example with `MatSetValues()` unless one first calls `MatSetUnfactored()`.
This is really in-place only for dense matrices, the preferred approach is to use `MatGetFactor()`, `MatLUFactorSymbolic()`, and `MatLUFactorNumeric()`
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), [Matrix Factorization](sec_matfactor), `Mat`, `MatFactorType`, `MatLUFactorSymbolic()`, `MatLUFactorNumeric()`, `MatCholeskyFactor()`,
`MatGetOrdering()`, `MatSetUnfactored()`, `MatFactorInfo`, `MatGetFactor()`
## Level
developer
## Location
src/mat/interface/matrix.c
## Examples
src/ksp/ksp/tutorials/ex74.c
## Implementations
MatLUFactor_SeqAIJ in src/mat/impls/aij/seq/aijfact.c
MatLUFactor_SeqBAIJ in src/mat/impls/baij/seq/baijfact.c
MatLUFactor_SeqDense in src/mat/impls/dense/seq/dense.c
MatLUFactor_Elemental in src/mat/impls/elemental/matelem.cxx
MatLUFactor_ScaLAPACK in src/mat/impls/scalapack/matscalapack.c
---
[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)