petsc-3.14.6 2021-03-30
Report Typos and Errors

TSComputeIJacobianDefaultColor

Computes the Jacobian using finite differences and coloring to exploit matrix sparsity.

Synopsis

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

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
ctx - an optional user context

Output Parameters

J - Jacobian matrix (not altered in this routine)
B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)

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.

This will first try to get the coloring from the DM. If the DM type has no coloring routine, then it will try to get the coloring from the matrix. This requires that the matrix have nonzero entries precomputed.

See Also

TSSetIJacobian(), MatFDColoringCreate(), MatFDColoringSetFunction()

Level

intermediate

Location

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