TSRKRegister#

register an TSRK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation

Synopsis#

#include "petscts.h"   
PetscErrorCode TSRKRegister(TSRKType name, PetscInt order, PetscInt s, const PetscReal A[], const PetscReal b[], const PetscReal c[], const PetscReal bembed[], PetscInt p, const PetscReal binterp[])

Not Collective, but the same schemes should be registered on all processes on which they will be used

Input Parameters#

  • name - identifier for method

  • order - approximation order of method

  • 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; NULL to use last row of A)

  • c - abscissa (dimension s; NULL to use row sums of A)

  • 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; NULL to reuse b with p=1)

Note#

Several TSRK methods are provided, this function is only needed to create new methods.

See Also#

TS: Scalable ODE and DAE Solvers, TSRK

Level#

advanced

Location#

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


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