TSComputeIJacobian#

Evaluates the Jacobian of the DAE

Synopsis#

#include "petscts.h"  
PetscErrorCode TSComputeIJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal shift, Mat A, Mat B, PetscBool imex)

Collective

Input

Input Parameters#

  • ts - the TS context

  • t - current timestep

  • U - state vector

  • Udot - time derivative of state vector

  • shift - shift to apply, see note below

  • imex - flag indicates if the method is TSIMEX so that the RHSJacobian should be kept separate

Output Parameters#

  • A - Jacobian matrix

  • B - matrix from which the preconditioner is constructed; often the same as A

Notes#

If F(t,U,Udot)=0 is the DAE, the required Jacobian is

   dF/dU + shift*dF/dUdot

Most users should not need to explicitly call this routine, as it is used internally within the nonlinear solvers.

See Also#

TS: Scalable ODE and DAE Solvers, TS, TSSetIJacobian()

Level#

developer

Location#

src/ts/interface/ts.c


Edit on GitLab

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