Manage coefficients of a functional (extended) Butcher tableau.
More...
|
| FunctionalButcherTableau ()=default |
| No memory allocation.
|
|
| FunctionalButcherTableau (unsigned s, unsigned order, const function_type *a, const function_type *b, const real_type *c) |
| Construct a classic non-embedded tableau.
|
|
| FunctionalButcherTableau (unsigned s, unsigned embedded_order, unsigned order, const function_type *a, const function_type *b, const function_type *bt, const real_type *c) |
| Construct an embedded tableau.
|
|
function_type | a (unsigned i, unsigned j) const |
| Read the a_ij coefficients.
|
|
real_type | c (unsigned i) const |
| Read the c_i coefficients.
|
|
function_type | b (unsigned j) const |
| Read the b_j coefficients.
|
|
function_type | bt (unsigned j) const |
| Read the embedded bt_j coefficients.
|
|
unsigned | num_stages () const |
| The number of stages s.
|
|
unsigned | order () const |
| global order of accuracy for the method represented by b
|
|
unsigned | embedded_order () const |
| global order of accuracy for the embedded method represented by bt
|
|
bool | isEmbedded () const |
| True if the method has an embedding.
|
|
bool | isImplicit () const |
| True if an element on or above the diagonal in a is non-zero.
|
|
template<class real_type>
struct dg::mat::FunctionalButcherTableau< real_type >
Manage coefficients of a functional (extended) Butcher tableau.
The goal of this class is to represent a Butcher tableau for the use in Exponential Runge Kutta type time integrators. The coefficients of the tableau are easily constructible and accessible. Furthermore, we provide utilities like the number of stages, whether the tableau is embedded or not and the order of the method.
Currently available are
We follow the naming convention of the ARKode library https://sundials.readthedocs.io/en/latest/arkode/Butcher_link.html (They also provide nice stability plots for their methods) as NAME-S-P-Q or NAME-S-Q, where
- NAME is the author or name of the method
- S is the number of stages in the method
- P is the global order of the embedding
- Q is the global order of the method
- Note
- In exponential Rosenbrock type schemes it is assumed that \( A\) (the matrix) is the Jacobian of the system. If it is not, then the order conditions are different and the order and embedded orders are not what is indicated in our names.
- Template Parameters
-
real_type | type of the coefficients |
- See also
- dg::mat::ExpRungeKutta
◆ function_type
template<class real_type >
◆ value_type
template<class real_type >
◆ FunctionalButcherTableau() [1/3]
template<class real_type >
◆ FunctionalButcherTableau() [2/3]
template<class real_type >
Construct a classic non-embedded tableau.
- Parameters
-
s | number of stages |
order | (global) order of the resulting method |
a | pointer to s*s real functions interpreted as a_{ij}=a[i*s+j] |
b | pointer to s real functions interpreted as b_j=b[j] |
c | pointer to s real numbers interpreted as c_i=c[i] |
- Note
- This constructor initializes the embedded coefficients bt=b which makes the embedded method equal to the actual method
◆ FunctionalButcherTableau() [3/3]
template<class real_type >
Construct an embedded tableau.
- Parameters
-
s | number of stages |
embedded_order | (global) order of the embedded method (corresponding to bt ) |
order | (global) order of the method (corresponding to b ) |
a | pointer to s*s real functions interpreted as a_{ij}=a[i*s+j] |
b | pointer to s real functions interpreted as b_j=b[j] |
bt | pointer to s real functions interpreted as bt_j=bt[j] |
c | pointer to s real numbers interpreted as c_i=c[i] |
◆ a()
template<class real_type >
Read the a_ij coefficients.
- Parameters
-
i | row number 0<=i<s, i>=s results in undefined behaviour |
j | col number 0<=j<s, j>=s results in undefined behaviour |
- Returns
- a_ij
- Note
- The returned function is either strictly positive or zero everywhere, so a test
a(i,j)(0) == 0
is valid
◆ b()
template<class real_type >
Read the b_j coefficients.
- Parameters
-
j | col number 0<=j<s, j>=s results in undefined behaviour |
- Returns
- b_j
- Note
- The returned function is either strictly positive or zero everywhere, so a test
b(j)(0) == 0
is valid
◆ bt()
template<class real_type >
Read the embedded bt_j coefficients.
- Parameters
-
j | col number 0<=j<s, j>=s results in undefined behaviour |
- Returns
- bt_j
- Note
- The returned function is either strictly positive or zero everywhere, so a test
bt(j)(0) == 0
is valid
◆ c()
template<class real_type >
Read the c_i coefficients.
- Parameters
-
i | row number 0<=i<s, i>=s results in undefined behaviour |
- Returns
- c_i
◆ embedded_order()
template<class real_type >
global order of accuracy for the embedded method represented by bt
◆ isEmbedded()
template<class real_type >
True if the method has an embedding.
◆ isImplicit()
template<class real_type >
True if an element on or above the diagonal in a is non-zero.
◆ num_stages()
template<class real_type >
◆ order()
template<class real_type >
global order of accuracy for the method represented by b
The documentation for this struct was generated from the following file: