petsc-3.14.6 2021-03-30
Report Typos and Errors

DMTSSetRHSJacobian

set TS Jacobian evaluation function

Synopsis

#include "petscts.h" 
PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx)
Not Collective

Input Argument

dm - DM to be used with TS
func - Jacobian evaluation routine
ctx - context for residual evaluation

Calling sequence of func

    PetscErrorCode func(TS ts,PetscReal t,Vec u,Mat A,Mat B,void *ctx);

t - current timestep
u - input vector
Amat - (approximate) Jacobian matrix
Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat)
ctx - [optional] user-defined context for matrix evaluation routine

Note

TSSetJacobian() is normally used, but it calls this function internally because the user context is actually associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.

See Also

DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetRHSJacobian()

Level

advanced

Location

src/ts/utils/dmts.c
Index of all TS routines
Table of Contents for all manual pages
Index of all manual pages