:orphan: # DMTSSetIJacobian set `TS` Jacobian evaluation function ## Synopsis ``` #include "petscts.h" PetscErrorCode DMTSSetIJacobian(DM dm, TSIJacobian func, void *ctx) ``` Not Collective ## Input Parameters - ***dm -*** `DM` to be used with `TS` - ***func -*** Jacobian evaluation routine - ***ctx -*** context for residual evaluation ## Calling sequence of `f` ```none PetscErrorCode f(TS ts, PetscReal t, Vec U, Vec U_t, PetscReal a, Mat Amat, Mat Pmat, void *ctx); ``` - ***ts -*** the `TS` context obtained from `TSCreate()` - ***t -*** time at step/stage being solved - ***U -*** state vector - ***U_t -*** time derivative of state vector - ***a -*** shift - ***Amat -*** (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t - ***Pmat -*** matrix used for constructing preconditioner, usually the same as Amat - ***ctx -*** [optional] user-defined context for matrix evaluation routine ## Note `TSSetJacobian()` is normally used, but it calls this function internally because the user context is actually associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or not. If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian. ## See Also [](ch_ts), `TS`, `DM`, `DMTSSetContext()`, `TSSetRHSFunction()`, `DMTSGetJacobian()`, `TSSetIJacobian()`, `TSSetIFunction()` ## Level advanced ## Location src/ts/utils/dmts.c --- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/ts/utils/dmts.c) [Index of all TS routines](index.md) [Table of Contents for all manual pages](/manualpages/index.md) [Index of all manual pages](/manualpages/singleindex.md)