TSGLEERegister#

register a new TSGLEE scheme by providing the entries in the Butcher tableau

Synopsis#

#include "petscts.h"   
PetscErrorCode TSGLEERegister(TSGLEEType name, PetscInt order, PetscInt s, PetscInt r, PetscReal gamma, const PetscReal A[], const PetscReal B[], const PetscReal U[], const PetscReal V[], const PetscReal S[], const PetscReal F[], const PetscReal c[], const PetscReal Fembed[], const PetscReal Ferror[], const PetscReal Serror[], PetscInt pinterp, 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 - order of method

  • s - number of stages

  • r - number of steps

  • gamma - LTE ratio

  • A - stage coefficients (dimension s*s, row-major)

  • B - step completion coefficients (dimension r*s, row-major)

  • U - method coefficients (dimension s*r, row-major)

  • V - method coefficients (dimension r*r, row-major)

  • S - starting coefficients

  • F - finishing coefficients

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

  • Fembed - step completion coefficients for embedded method

  • Ferror - error computation coefficients

  • Serror - error initialization coefficients

  • pinterp - order of interpolation (0 if unavailable)

  • binterp - array of interpolation coefficients (NULL if unavailable)

Note#

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

See Also#

TS: Scalable ODE and DAE Solvers, TSGLEE

Level#

advanced

Location#

src/ts/impls/glee/glee.c


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