TSDIRKRegister#
register a TSDIRK
scheme by providing the entries in its Butcher tableau and, optionally, embedded approximations and interpolation
Synopsis#
#include "petscts.h"
PetscErrorCode TSDIRKRegister(TSDIRKType name, PetscInt order, PetscInt s, const PetscReal At[], const PetscReal bt[], const PetscReal ct[], const PetscReal bembedt[], PetscInt pinterp, const PetscReal binterpt[])
Logically Collective.
Input Parameters#
name - identifier for method
order - approximation order of method
s - number of stages, this is the dimension of the matrices below
At - Butcher table of stage coefficients (dimension
s
*s
, row-major order)bt - Butcher table for completing the step (dimension
s
; passNULL
to use the last row ofAt
)ct - Abscissa of each stage (dimension s, NULL to use row sums of At)
bembedt - Stiff part of completion table for embedded method (dimension s;
NULL
if not available)pinterp - Order of the interpolation scheme, equal to the number of columns of
binterpt
andbinterp
binterpt - Coefficients of the interpolation formula (dimension s*pinterp)
Note#
Several TSDIRK
methods are provided, the use of this function is only needed to create new methods.
See Also#
Level#
advanced
Location#
src/ts/impls/arkimex/arkimex.c
Index of all TS routines
Table of Contents for all manual pages
Index of all manual pages