TSEvaluateWLTE#

Evaluate the weighted local truncation error norm at the end of a time step with a given order of accuracy.

Synopsis#

#include "petscts.h"  
PetscErrorCode TSEvaluateWLTE(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte)

Collective

Input Parameters#

Input/Output Parameter#

  • order - optional, desired order for the error evaluation or PETSC_DECIDE; on output, the actual order of the error evaluation

Output Parameter#

  • wlte - the weighted local truncation error norm

Note#

If the timestepper cannot evaluate the error in a particular step (eg. in the first step or restart steps after event handling), this routine returns wlte=-1.0 .

See Also#

TS: Scalable ODE and DAE Solvers, TS, TSStep(), TSAdapt, TSErrorWeightedNorm()

Level#

advanced

Location#

src/ts/interface/ts.c

Implementations#

TSEvaluateWLTE_BDF in src/ts/impls/bdf/bdf.c
TSEvaluateWLTE_Alpha in src/ts/impls/implicit/alpha/alpha1.c
TSEvaluateWLTE_Alpha in src/ts/impls/implicit/alpha/alpha2.c
TSEvaluateWLTE_Theta in src/ts/impls/implicit/theta/theta.c


Edit on GitLab

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