2.6. SPRK Method Table Structure
To store the pair of Butcher tables defining a SPRK method of order \(q\)
ARKODE provides the ARKodeSPRKTable
type and several related utility
routines. We use the following notation
where \(B\) and \(\hat{B}\) contain the coefficients for the explicit
and diagonally implicit tables, respectively. We use a compact storage of these
coefficients in terms of two arrays, one for \(a\) values and one for
\(\hat{a}\) values. The abscissae (only relevant for non-autonomous
problems) are computed dynamically as \(c_j = \sum_{i=1}^j a_i\) and
\(\hat{c}_j = \sum_{i=1}^j \hat{a}_i\), respectively
[44, 77]. The ARKodeSPRKTable
type is a pointer to
the ARKodeSPRKTableMem
structure:
-
typedef ARKodeSPRKTableMem *ARKodeSPRKTable
-
type ARKodeSPRKTableMem
Structure representing the SPRK method that holds the method coefficients.
-
int q
The method order of accuracy.
-
int stages
The number of stages.
-
sunrealtype *a
Array of coefficients that generate the explicit Butcher table.
a[i]
is the coefficient appearing in column i+1.
-
sunrealtype *ahat
Array of coefficients that generate the diagonally-implicit Butcher table.
ahat[i]
is the coefficient appearing in column i.
-
int q
2.6.1. ARKodeSPRKTable functions
Function name |
Description |
---|---|
Allocate an empty table |
|
Load SPRK method using an identifier |
|
Load SPRK method using a string version of the identifier |
|
Create a new table |
|
Create a copy of a table |
|
Get the table real and integer workspace size |
|
Deallocate a table |
-
ARKodeSPRKTable ARKodeSPRKTable_Create(int stages, int q, const sunrealtype *a, const sunrealtype *ahat)
Creates and allocates an
ARKodeSPRKTable
with the specified number of stages and the coefficients provided.- Parameters:
stages – The number of stages.
q – The order of the method.
a – An array of the coefficients for the
a
table.ahat – An array of the coefficients for the
ahat
table.
- Returns:
ARKodeSPRKTable
for the loaded method.
-
ARKodeSPRKTable ARKodeSPRKTable_Alloc(int stages)
Allocate memory for an
ARKodeSPRKTable
with the specified number of stages.- Parameters:
stages – The number of stages.
- Returns:
ARKodeSPRKTable
for the loaded method.
-
ARKodeSPRKTable ARKodeSPRKTable_Load(ARKODE_SPRKMethodID id)
Load the
ARKodeSPRKTable
for the specified method ID.- Parameters:
id – The ID of the SPRK method, see Symplectic Partitioned Butcher tables.
- Returns:
ARKodeSPRKTable
for the loaded method.
-
ARKodeSPRKTable ARKodeSPRKTable_LoadByName(const char *method)
Load the
ARKodeSPRKTable
for the specified method name.- Parameters:
method – The name of the SPRK method, see Symplectic Partitioned Butcher tables.
- Returns:
ARKodeSPRKTable
for the loaded method.
-
ARKodeSPRKTable ARKodeSPRKTable_Copy(ARKodeSPRKTable sprk_table)
Create a copy of the
ARKodeSPRKTable
.- Parameters:
sprk_table – The
ARKodeSPRKTable
to copy.
- Returns:
Pointer to the copied
ARKodeSPRKTable
.
-
void ARKodeSPRKTable_Write(ARKodeSPRKTable sprk_table, FILE *outfile)
Write the ARKodeSPRKTable out to the file.
- Parameters:
sprk_table – The
ARKodeSPRKTable
to write.outfile – The FILE that will be written to.
-
void ARKodeSPRKTable_Space(ARKodeSPRKTable sprk_table, sunindextype *liw, sunindextype *lrw)
Get the workspace sizes required for the
ARKodeSPRKTable
.- Parameters:
sprk_table – The
ARKodeSPRKTable
.liw – Pointer to store the integer workspace size.
lrw – Pointer to store the real workspace size.
-
void ARKodeSPRKTable_Free(ARKodeSPRKTable sprk_table)
Free the memory allocated for the
ARKodeSPRKTable
.- Parameters:
sprk_table – The
ARKodeSPRKTable
to free.
-
int ARKodeSPRKTable_ToButcher(ARKodeSPRKTable sprk_table, ARKodeButcherTable *a_ptr, ARKodeButcherTable *b_ptr)
Convert the
ARKodeSPRKTable
to the Butcher table representation.- Parameters:
sprk_table – The
ARKodeSPRKTable
.a_ptr – Pointer to store the explicit Butcher table.
b_ptr – Pointer to store the diagonally-implicit Butcher table.