:orphan:
# MatCreateMFFD
Creates a matrix-free matrix of type `MATMFFD`. See also `MatCreateSNESMF()` 
## Synopsis
```
#include "petscmat.h"   
PetscErrorCode MatCreateMFFD(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, Mat *J)
```
Collective


## Input Parameters

- ***comm -*** MPI communicator
- ***m -*** number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given)
This value should be the same as the local size used in creating the
y vector for the matrix-vector product y = Ax.
- ***n -*** This value should be the same as the local size used in creating the
x vector for the matrix-vector product y = Ax. (or `PETSC_DECIDE` to have
calculated if `N` is given) For square matrices `n` is almost always `m`.
- ***M -*** number of global rows (or `PETSC_DETERMINE` to have calculated if `m` is given)
- ***N -*** number of global columns (or `PETSC_DETERMINE` to have calculated if `n` is given)



## Output Parameter

- ***J -*** the matrix-free matrix



## Options Database Keys

- ***-mat_mffd_type -*** wp or ds (see `MATMFFD_WP` or `MATMFFD_DS`)
- ***-mat_mffd_err -*** square root of estimated relative error in function evaluation
- ***-mat_mffd_period -*** how often h is recomputed, defaults to 1, every time
- ***-mat_mffd_check_positivity -*** possibly decrease `h` until U + h*a has only positive values
- ***-mat_mffd_umin <umin> -*** Sets umin (for default PETSc routine that computes h only)
- ***-mat_mffd_complex -*** use the Lyness trick with complex numbers to compute the matrix-vector product instead of differencing
(requires real valued functions but that PETSc be configured for complex numbers)





## Notes
The matrix-free matrix context merely contains the function pointers
and work space for performing finite difference approximations of
Jacobian-vector products, F'(u)*a,

The default code uses the following approach to compute h

```none
     F'(u)*a = [F(u+h*a) - F(u)]/h where
     h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
       = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
 where
     error_rel = square root of relative error in function evaluation
     umin = minimum iterate parameter
```


You can call `SNESSetJacobian()` with `MatMFFDComputeJacobian()` if you are using matrix and not a different
preconditioner matrix

The user can set the error_rel via `MatMFFDSetFunctionError()` and
umin via `MatMFFDDSSetUmin()`.

The user should call `MatDestroy()` when finished with the matrix-free
matrix context.


## See Also
 [](ch_matrices), `Mat`, `MATMFFD`, `MatDestroy()`, `MatMFFDSetFunctionError()`, `MatMFFDDSSetUmin()`, `MatMFFDSetFunction()`
`MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `MatCreateSNESMF()`,
`MatMFFDGetH()`, `MatMFFDRegister()`, `MatMFFDComputeJacobian()`

## Level
advanced

## Location
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/mat/impls/mffd/mffd.c.html#MatCreateMFFD">src/mat/impls/mffd/mffd.c</A>

## Examples
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/snes/tutorials/ex22.c.html">src/snes/tutorials/ex22.c</A><BR>

## Implementations

<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/mat/impls/mffd/mffddef.c.html#MatCreateMFFD_DS">MatCreateMFFD_DS in src/mat/impls/mffd/mffddef.c</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/mat/impls/mffd/wp.c.html#MatCreateMFFD_WP">MatCreateMFFD_WP in src/mat/impls/mffd/wp.c</A><BR>


---
[Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/mat/impls/mffd/mffd.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)  
