petsc-3.7.7 2017-09-25
Report Typos and Errors

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 on TS

Input Arguments

ts - time stepping context
wnormtype - norm type, either NORM_2 or NORM_INFINITY
order - optional, desired order for the error evaluation or PETSC_DECIDE

Output Arguments

order - optional, the actual order of the error evaluation
wlte - the weighted local truncation error norm

Notes

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

TSStep(), TSAdapt, TSErrorWeightedNorm()

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