#include "petscts.h" PetscErrorCode TSARKIMEXRegister(TSARKIMEXType name,PetscInt order,PetscInt s, const PetscReal At[],const PetscReal bt[],const PetscReal ct[], const PetscReal A[],const PetscReal b[],const PetscReal c[], const PetscReal bembedt[],const PetscReal bembed[], PetscInt pinterp,const PetscReal binterpt[],const PetscReal binterp[])Not Collective, but the same schemes should be registered on all processes on which they will be used
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 for stiff part (dimension s*s, row-major) | |
bt | - Butcher table for completing the stiff part of the step (dimension s; NULL to use the last row of At) | |
ct | - Abscissa of each stiff stage (dimension s, NULL to use row sums of At) | |
A | - Non-stiff stage coefficients (dimension s*s, row-major) | |
b | - Non-stiff step completion table (dimension s; NULL to use last row of At) | |
c | - Non-stiff abscissa (dimension s; NULL to use row sums of A) | |
bembedt | - Stiff part of completion table for embedded method (dimension s; NULL if not available) | |
bembed | - Non-stiff part of completion table for embedded method (dimension s; NULL to use bembedt if provided) | |
pinterp | - Order of the interpolation scheme, equal to the number of columns of binterpt and binterp | |
binterpt | - Coefficients of the interpolation formula for the stiff part (dimension s*pinterp) | |
binterp | - Coefficients of the interpolation formula for the non-stiff part (dimension s*pinterp; NULL to reuse binterpt) |
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