|
| FunctionalButcherTableau ()=default |
| No memory allocation. More...
|
|
| FunctionalButcherTableau (unsigned s, unsigned order, const function_type *a, const function_type *b, const real_type *c) |
| Construct a classic non-embedded tableau. More...
|
|
| 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. More...
|
|
function_type | a (unsigned i, unsigned j) const |
| Read the a_ij coefficients. More...
|
|
real_type | c (unsigned i) const |
| Read the c_i coefficients. More...
|
|
function_type | b (unsigned j) const |
| Read the b_j coefficients. More...
|
|
function_type | bt (unsigned j) const |
| Read the embedded bt_j coefficients. More...
|
|
unsigned | num_stages () const |
| The number of stages s. More...
|
|
unsigned | order () const |
| global order of accuracy for the method represented by b More...
|
|
unsigned | embedded_order () const |
| global order of accuracy for the embedded method represented by bt More...
|
|
bool | isEmbedded () const |
| True if the method has an embedding. More...
|
|
bool | isImplicit () const |
| True if an element on or above the diagonal in a is non-zero. More...
|
|
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 http://runge.math.smu.edu/arkode_dev/doc/guide/build/html/Butcher.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