DMTSSetRHSJacobian#

set TS Jacobian evaluation function

Synopsis#

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

Not Collective

Input Parameters#

  • 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);
  • ts - the TS context obtained from TSCreate()

  • 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.

Developer Note#

If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.

See Also#

TS: Scalable ODE and DAE Solvers, DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetRHSJacobian()

Level#

advanced

Location#

src/ts/utils/dmts.c


Edit on GitLab

Index of all TS routines
Table of Contents for all manual pages
Index of all manual pages