petsc-3.7.1 2016-05-15
Report Typos and Errors

TSSetIJacobian

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

Synopsis

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

Input Parameters

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

Calling sequence of f

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

t - time at step/stage being solved Many br
U - state vector Many br
U_t - time derivative of state vector Many br
a - shift Many br
Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t Many br
Pmat - matrix used for constructing preconditioner, usually the same as Amat Many br
ctx - [optional] user-defined context for matrix evaluation routine Many br

Notes

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

If you know the operator Amat has a null space you can use MatSetNullSpace() and MatSetTransposeNullSpace() to supply the null Many brspace to Amat and the KSP solvers will automatically use that null space as needed during the solution process. Many br

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

You must set all the diagonal entries of the matrices, if they are zero you must still set them with a zero value Many br

The TS solver may modify the nonzero structure and the entries of the matrices Amat and Pmat between the calls to f() Many brYou should not assume the values are the same in the next call to f() as you set them in the previous call. Many br

Many br

Keywords

TS, timestep, DAE, Jacobian

See Also

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

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/ex19.c.html
src/ts/examples/tutorials/ex20.c.html
src/ts/examples/tutorials/ex22.c.html