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, TSIJacobianFn *f, void *ctx)
Logically Collective
Input Parameters#
Amat - (approximate) matrix to store Jacobian entries computed by
f
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
)
Notes#
The matrices Amat
and Pmat
are exactly the matrices that are used by SNES
for the nonlinear solve.
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 + adF/dU_t you provide turns out to be the Jacobian of F(t,U,W+aU) where F(t,U,U_t) = 0 is the DAE to be solved. The time integrator internally approximates U_t by W+aU 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 = -aU(previous timestep) so W + aU = 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.
In case TSSetRHSJacobian()
is also used in conjunction with a fully-implicit solver,
multilevel linear solvers, e.g. PCMG
, will likely not work due to the way TS
handles rhs matrices.
See Also#
TS: Scalable ODE and DAE Solvers, TS
, TSIJacobianFn
, TSSetIFunction()
, TSSetRHSJacobian()
,
SNESComputeJacobianDefaultColor()
, SNESComputeJacobianDefault()
, TSSetRHSFunction()
Level#
beginner
Location#
Examples#
src/ts/tutorials/ex31.c
src/ts/tutorials/ex20adj.c
src/ts/tutorials/ex20.c
src/ts/tutorials/ex20fwd.c
src/ts/utils/dmplexlandau/tutorials/ex1f90.F90
src/ts/tutorials/ex32.c
src/ts/tutorials/ex24.c
src/ts/utils/dmplexlandau/tutorials/ex2.c
src/ts/utils/dmplexlandau/tutorials/ex1.c
src/ts/tutorials/ex36.c
Index of all TS routines
Table of Contents for all manual pages
Index of all manual pages