petsc-3.7.7 2017-09-25
Report Typos and Errors

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 on TS and Vec

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 IMEX so that the RHSJacobian should be kept separate

Output Parameters

A - Jacobian matrix
B - optional preconditioning matrix
flag - flag indicating matrix structure

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.

Keywords

TS, compute, Jacobian, matrix

See Also

TSSetIJacobian()

Level:developer
Location:
src/ts/interface/ts.c
Index of all TS routines
Table of Contents for all manual pages
Index of all manual pages