petsc-3.4.5 2014-06-29

TSSetIJacobian

Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function you provided with TSSetIFunction().

Synopsis

#include "petscts.h"  
PetscErrorCode  TSSetIJacobian(TS ts,Mat Amat,Mat Pmat,TSIJacobian f,void *ctx)
Logically Collective on TS

Input Parameters

ts - the TS context obtained from TSCreate()
Amat - (approximate) Jacobian matrix
Pmat - matrix used to compute preconditioner (usually the same as Amat)
f - the Jacobian evaluation routine
ctx - user-defined context for private data for the Jacobian evaluation routine (may be NULL)

Calling sequence of f

 f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat *Amat,Mat *Pmat,MatStructure *flag,void *ctx);

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
flag - flag indicating information about the preconditioner matrix structure (same as flag in KSPSetOperators())
ctx - [optional] user-defined context for matrix evaluation routine

Notes

The matrices Amat and Pmat are exactly the matrices that are used by SNES for the nonlinear solve.

The matrix dF/dU + a*dF/dU_t you provide turns out to be the Jacobian of F(t,U,W+a*U) where F(t,U,U_t) = 0 is the DAE to be solved. The time integrator internally approximates U_t by W+a*U where the positive "shift" a and vector W depend on the integration method, step size, and past states. For example with the backward Euler method a = 1/dt and W = -a*U(previous timestep) so W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt

Keywords

TS, timestep, DAE, Jacobian

See Also

TSSetIFunction(), TSSetRHSJacobian(), SNESComputeJacobianDefaultColor(), SNESComputeJacobianDefault()

Level:beginner
Location:
src/ts/interface/ts.c
Index of all TS routines
Table of Contents for all manual pages
Index of all manual pages

Examples

src/ts/examples/tutorials/ex6.c.html
src/ts/examples/tutorials/ex8.c.html
src/ts/examples/tutorials/ex9.c.html
src/ts/examples/tutorials/ex10.c.html
src/ts/examples/tutorials/ex14.c.html
src/ts/examples/tutorials/ex15.c.html
src/ts/examples/tutorials/ex16.c.html
src/ts/examples/tutorials/ex17.c.html
src/ts/examples/tutorials/ex22.c.html
src/ts/examples/tutorials/ex23.c.html
src/ts/examples/tutorials/ex24.c.html