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

TSRKGetTableau

Get info on the RK tableau

Synopsis

#include "petscts.h"   
PetscErrorCode TSRKGetTableau(TS ts, PetscInt *s, const PetscReal **A, const PetscReal **b, const PetscReal **c, const PetscReal **bembed,
                                     PetscInt *p, const PetscReal **binterp, PetscBool *FSAL)
Not Collective

Input Parameters

ts - timestepping context

Output Parameters

s - number of stages, this is the dimension of the matrices below
A - stage coefficients (dimension s*s, row-major)
b - step completion table (dimension s)
c - abscissa (dimension s)
bembed - completion table for embedded method (dimension s; NULL if not available)
p - Order of the interpolation scheme, equal to the number of columns of binterp
binterp - Coefficients of the interpolation formula (dimension s*p)
FSAL - wheather or not the scheme has the First Same As Last property

See Also

TSRK

Level

developer

Location

src/ts/impls/explicit/rk/rk.c

Implementations

src/ts/impls/explicit/rk/rk.c:463:PetscErrorCode TSRKGetTableau_RK(TS ts, PetscInt *s, const PetscReal **A, const PetscReal **b, const PetscReal **c, const PetscReal **bembed,
Index of all TS routines
Table of Contents for all manual pages
Index of all manual pages