petsc-3.7.1 2016-05-15
Report Typos and Errors

TSRKRegister

register an RK 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 pinterp,const PetscReal binterp[])
Not Collective, but the same schemes should be registered on all processes on which they will be used Many br

Input Parameters

name - identifier for method Many br
order - approximation order of method Many br
s - number of stages, this is the dimension of the matrices below Many br
A - stage coefficients (dimension s*s, row-major) Many br
b - step completion table (dimension s; NULL to use last row of A) Many br
c - abscissa (dimension s; NULL to use row sums of A) Many br
bembed - completion table for embedded method (dimension s; NULL if not available) Many br
pinterp - Order of the interpolation scheme, equal to the number of columns of binterp Many br
binterp - Coefficients of the interpolation formula (dimension s*pinterp; NULL to reuse binterpt) Many br

Notes

Several RK methods are provided, this function is only needed to create new methods. Many br

Many br

Keywords

TS, register

See Also

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