TSSetI2Jacobian#

Set the function to compute the matrix dF/dU + vdF/dU_t + adF/dU_tt where F(t,U,U_t,U_tt) is the function you provided with TSSetI2Function().

Synopsis#

#include "petscts.h"  
PetscErrorCode TSSetI2Jacobian(TS ts, Mat J, Mat P, TSI2JacobianFn *jac, void *ctx)

Logically Collective

Input Parameters#

  • ts - the TS context obtained from TSCreate()

  • J - matrix to hold the Jacobian values

  • P - matrix for constructing the preconditioner (may be same as J)

  • jac - the Jacobian evaluation routine, see TSI2JacobianFn for the calling sequence

  • ctx - user-defined context for private data for the Jacobian evaluation routine (may be NULL)

Notes#

The matrices J and P are exactly the matrices that are used by SNES for the nonlinear solve.

The matrix dF/dU + vdF/dU_t + adF/dU_tt you provide turns out to be the Jacobian of G(U) = F(t,U,W+vU,W’+aU) where F(t,U,U_t,U_tt) = 0 is the DAE to be solved. The time integrator internally approximates U_t by W+vU and U_tt by W’+aU where the positive “shift” parameters ‘v’ and ‘a’ and vectors W, W’ depend on the integration method, step size, and past states.

See Also#

TS: Scalable ODE and DAE Solvers, TS, TSI2JacobianFn, TSSetI2Function(), TSGetI2Jacobian()

Level#

beginner

Location#

src/ts/interface/ts.c

Examples#

src/ts/tutorials/ex44.c
src/ts/tutorials/ex43.c


Index of all TS routines
Table of Contents for all manual pages
Index of all manual pages