# TSComputeI2Jacobian
Evaluates the Jacobian of the DAE 
## Synopsis
```
#include "petscts.h"  
PetscErrorCode TSComputeI2Jacobian(TS ts, PetscReal t, Vec U, Vec V, Vec A, PetscReal shiftV, PetscReal shiftA, Mat J, Mat P)
```
Collective


## Input Parameters

- ***ts -*** the `TS` context
- ***t -*** current timestep
- ***U -*** state vector
- ***V -*** time derivative of state vector
- ***A -*** second time derivative of state vector
- ***shiftV -*** shift to apply, see note below
- ***shiftA -*** shift to apply, see note below



## Output Parameters

- ***J -*** Jacobian matrix
- ***P -*** optional preconditioning matrix





## Notes
If F(t,U,V,A)=0 is the DAE, the required Jacobian is

dF/dU + shiftV*dF/dV + shiftA*dF/dA

Most users should not need to explicitly call this routine, as it
is used internally within the nonlinear solvers.


## See Also
 [](chapter_ts), `TS`, `TSSetI2Jacobian()`

## Level
developer

## Location
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/interface/ts.c.html#TSComputeI2Jacobian">src/ts/interface/ts.c</A>


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


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