petsc-3.10.5 2019-03-28
Report Typos and Errors

TSComputeIJacobianConstant

Reuses a time-independent for a semi-implicit DAE or ODE

Synopsis

#include "petscts.h"  
PetscErrorCode TSComputeIJacobianConstant(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat A,Mat B,void *ctx)
Collective on TS

Input Arguments

ts - time stepping context
t - time at which to evaluate
U - state at which to evaluate
Udot - time derivative of state vector
shift - shift to apply
ctx - context

Output Arguments

A - pointer to operator
B - pointer to preconditioning matrix
flg - matrix structure flag

Notes

This function is intended to be passed to TSSetIJacobian() to evaluate the Jacobian for linear time-independent problems.

It is only appropriate for problems of the form

    M Udot = F(U,t)

where M is constant and F is non-stiff. The user must pass M to TSSetIJacobian(). The current implementation only works with IMEX time integration methods such as TSROSW and TSARKIMEX, since there is no support for de-constructing an implicit operator of the form

   shift*M + J

where J is the Jacobian of -F(U). Support may be added in a future version of PETSc, but for now, the user must store a copy of M or reassemble it when requested.

See Also

TSSetIFunction(), TSSetIJacobian(), TSComputeIFunctionLinear()

Level

advanced

Location

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