#include "petscts.h" PetscErrorCode TSSetI2Jacobian(TS ts,Mat J,Mat P,TSI2Jacobian jac,void *ctx)Logically Collective on TS
ts | - the TS context obtained from TSCreate() | |
J | - Jacobian matrix | |
P | - preconditioning matrix for J (may be same as J) | |
jac | - the Jacobian evaluation routine | |
ctx | - user-defined context for private data for the Jacobian evaluation routine (may be NULL) |
jac(TS ts,PetscReal t,Vec U,Vec U_t,Vec U_tt,PetscReal v,PetscReal a,Mat J,Mat P,void *ctx);
t | - time at step/stage being solved | |
U | - state vector | |
U_t | - time derivative of state vector | |
U_tt | - second time derivative of state vector | |
v | - shift for U_t | |
a | - shift for U_tt | |
J | - Jacobian of G(U) = F(t,U,W+v*U,W'+a*U), equivalent to dF/dU + v*dF/dU_t + a*dF/dU_tt | |
P | - preconditioning matrix for J, may be same as J | |
ctx | - [optional] user-defined context for matrix evaluation routine |
The matrix dF/dU + v*dF/dU_t + a*dF/dU_tt you provide turns out to be the Jacobian of G(U) = F(t,U,W+v*U,W'+a*U) where F(t,U,U_t,U_tt) = 0 is the DAE to be solved. The time integrator internally approximates U_t by W+v*U and U_tt by W'+a*U where the positive "shift" parameters 'v' and 'a' and vectors W, W' depend on the integration method, step size, and past states.
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