PETSc version 3.15.5
TSMPRKRegister
register a MPRK scheme by providing the entries in the Butcher tableau 
Synopsis
#include "petscts.h"   
PetscErrorCode TSMPRKRegister(TSMPRKType name,PetscInt order,
                              PetscInt sbase,PetscInt ratio1,PetscInt ratio2,
                              const PetscReal Asb[],const PetscReal bsb[],const PetscReal csb[],const PetscInt rsb[],
                              const PetscReal Amb[],const PetscReal bmb[],const PetscReal cmb[],const PetscInt rmb[],
                              const PetscReal Af[],const PetscReal bf[],const PetscReal cf[])
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 in the base methods | 
|  | ratio1 | - stepsize ratio at 1st level (e.g. slow/medium) | 
|  | ratio2 | - stepsize ratio at 2nd level (e.g. medium/fast) | 
|  | Af | - stage coefficients for fast components(dimension s*s, row-major) | 
|  | bf | - step completion table for fast components(dimension s) | 
|  | cf | - abscissa for fast components(dimension s) | 
|  | As | - stage coefficients for slow components(dimension s*s, row-major) | 
|  | bs | - step completion table for slow components(dimension s) | 
|  | cs | - abscissa for slow components(dimension s) | 
Notes
Several MPRK methods are provided, this function is only needed to create new methods.
See Also
 TSMPRK
Level
advanced
Location
src/ts/impls/multirate/mprk.c
Index of all TS routines
Table of Contents for all manual pages
Index of all manual pages