#include "petscts.h" PetscErrorCode TSSetRHSJacobian(TS ts,Mat A,Mat B,TSRHSJacobian f,void *ctx)Logically Collective on TS
ts | - the TS context obtained from TSCreate() | |
A | - Jacobian matrix | |
B | - preconditioner matrix (usually same as A) | |
f | - the Jacobian evaluation routine | |
ctx | - [optional] user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL) |
func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
t | - current timestep | |
u | - input vector | |
A | - matrix A, where U_t = A(t)u | |
B | - preconditioner matrix, usually the same as A | |
flag | - flag indicating information about the preconditioner matrix structure (same as flag in KSPSetOperators()) | |
ctx | - [optional] user-defined context for matrix evaluation routine |
The routine func() takes Mat * as the matrix arguments rather than Mat. This allows the matrix evaluation routine to replace A and/or B with a completely new matrix structure (not just different matrix elements) when appropriate, for instance, if the nonzero structure is changing throughout the global iterations.
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