petsc-3.9.4 2018-09-11
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
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
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
Index of all TS routines
Table of Contents for all manual pages
Index of all manual pages