petsc-3.7.3 2016-08-01
Report Typos and Errors

TSComputeI2Jacobian

Evaluates the Jacobian of the DAE

Synopsis

#include "petscts.h"  
PetscErrorCode TSComputeI2Jacobian(TS ts,PetscReal t,Vec U,Vec V,Vec A,PetscReal shiftV,PetscReal shiftA,Mat J,Mat P)
Collective on TS and Vec

Input Parameters

ts - the TS context
t - current timestep
U - state vector
V - time derivative of state vector
A - second time derivative of state vector
shiftV - shift to apply, see note below
shiftA - shift to apply, see note below

Output Parameters

J - Jacobian matrix
P - optional preconditioning matrix

Notes

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

dF/dU + shiftV*dF/dV + shiftA*dF/dA

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

TSSetI2Jacobian()

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