TSPseudoComputeFunction#

Compute nonlinear residual for pseudo-timestepping

Synopsis#

#include "petscts.h"   
PetscErrorCode TSPseudoComputeFunction(TS ts, Vec solution, Vec *residual, PetscReal *fnorm)

This computes the residual for \(\dot U = 0\), i.e. \(F(U, 0)\) for the IFunction.

Collective

Input Parameters#

  • ts - the timestep context

  • solution - the solution vector

Output Parameter#

  • residual - the nonlinear residual

  • fnorm - the norm of the nonlinear residual

Note#

TSPSEUDO records the nonlinear residual and the solution vector used to generate it. If given the same solution vector (as determined by the vectors PetscObjectState), this function will return those recorded values.

This would be used in a custom adaptive timestepping implementation that needs access to the residual, but reuses the calculation done by TSPSEUDO by default.

To correctly get the residual reuse behavior, solution must be the same Vec that returned by TSGetSolution().

See Also#

TS: Scalable ODE and DAE Solvers, TSPSEUDO

Level#

advanced

Location#

src/ts/impls/pseudo/posindep.c


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