#include "petscts.h" PetscErrorCode TSSetIJacobian(TS ts,Mat Amat,Mat Pmat,TSIJacobian f,void *ctx)Logically Collective on TS
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) |
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 | |
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 | |
ctx | - [optional] user-defined context for matrix evaluation routine |
If you know the operator Amat has a null space you can use MatSetNullSpace() and MatSetTransposeNullSpace() to supply the null space to Amat and the KSP solvers will automatically use that null space as needed during the solution process.
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
You must set all the diagonal entries of the matrices, if they are zero you must still set them with a zero value
The TS solver may modify the nonzero structure and the entries of the matrices Amat and Pmat between the calls to f() You should not assume the values are the same in the next call to f() as you set them in the previous call.