# MatMFFDComputeJacobian
Tells the matrix-free Jacobian object the new location at which Jacobian matrix vector products will be computed at, i.e. J(x) * a. The x is obtained from the `SNES` object (using `SNESGetSolution()`). 
## Synopsis
```
#include "petscsnes.h" 
#include "petscdm.h"   
PetscErrorCode MatMFFDComputeJacobian(SNES snes, Vec x, Mat jac, Mat B, void *dummy)
```
Collective


## Input Parameters

- ***snes -*** the nonlinear solver context
- ***x -*** the point at which the Jacobian vector products will be performed
- ***jac -*** the matrix-free Jacobian object of `MatType` `MATMFFD`, likely obtained with `MatCreateSNESMF()`
- ***B -*** either the same as jac or another matrix type (ignored)
- ***flag -*** not relevant for matrix-free form
- ***dummy -*** the user context (ignored)



## Option Database Key

- ***-snes_mf -*** use the matrix created with `MatSNESMFCreate()` to setup the Jacobian for each new solution in the Newton process





## Notes
If `MatMFFDSetBase()` is ever called on jac then this routine will NO longer get
the x from the `SNES` object and `MatMFFDSetBase()` must from that point on be used to
change the base vector x.

This can be passed into `SNESSetJacobian()` as the Jacobian evaluation function argument
when using a completely matrix-free solver,
that is the B matrix is also the same matrix operator. This is used when you select
-snes_mf but rarely used directly by users. (All this routine does is call `MatAssemblyBegin/End()` on
the `Mat` jac.)


## See Also
 `MatMFFDGetH()`, `MatCreateSNESMF()`, `MatMFFDSetBase()`, `MatCreateMFFD()`, `MATMFFD`,
`MatMFFDSetHHistory()`, `MatMFFDSetFunctionError()`, `MatCreateMFFD()`, `SNESSetJacobian()`

## Level
developer

## Location
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/snes/mf/snesmfj.c.html#MatMFFDComputeJacobian">src/snes/mf/snesmfj.c</A>


---
[Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/snes/mf/snesmfj.c)


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