Extension: Matrix functions
#include "dg/matrix/matrix.h"
dg::mat::FunctionalButcherTableau< real_type > Struct Template Reference

Manage coefficients of a functional (extended) Butcher tableau. More...

Public Types

using value_type = real_type
 
using function_type = std::function< value_type(value_type)>
 

Public Member Functions

 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...
 

Detailed Description

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
Name Identifier Description
Euler dg::mat::EXPLICIT_EULER_1_1 Hochbruck and Ostermann, Exponential Integrators, Acta Numerica (2010)
Midpoint-2-2 dg::mat::MIDPOINT_2_2 Hochbruck and Ostermann, Exponential Integrators, Acta Numerica (2010)
Runge-Kutta-4-4 dg::mat::CLASSIC_4_4 Cox and Matthews, J. Comput. Phys., 176 (2002)
Hochbruck-3-3-4 dg::mat::HOCHBRUCK_3_3_4 Hochbruck and Ostermann, Exponential Integrators, Acta Numerica (2010) (The exprb43 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_typetype of the coefficients
See also
dg::mat::ExpRungeKutta

Member Typedef Documentation

◆ function_type

template<class real_type >
using dg::mat::FunctionalButcherTableau< real_type >::function_type = std::function<value_type(value_type)>

◆ value_type

template<class real_type >
using dg::mat::FunctionalButcherTableau< real_type >::value_type = real_type

Constructor & Destructor Documentation

◆ FunctionalButcherTableau() [1/3]

template<class real_type >
dg::mat::FunctionalButcherTableau< real_type >::FunctionalButcherTableau ( )
default

No memory allocation.

◆ FunctionalButcherTableau() [2/3]

template<class real_type >
dg::mat::FunctionalButcherTableau< real_type >::FunctionalButcherTableau ( unsigned  s,
unsigned  order,
const function_type a,
const function_type b,
const real_type *  c 
)
inline

Construct a classic non-embedded tableau.

Parameters
snumber of stages
order(global) order of the resulting method
apointer to s*s real functions interpreted as a_{ij}=a[i*s+j]
bpointer to s real functions interpreted as b_j=b[j]
cpointer 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 >
dg::mat::FunctionalButcherTableau< real_type >::FunctionalButcherTableau ( unsigned  s,
unsigned  embedded_order,
unsigned  order,
const function_type a,
const function_type b,
const function_type bt,
const real_type *  c 
)
inline

Construct an embedded tableau.

Parameters
snumber of stages
embedded_order(global) order of the embedded method (corresponding to bt)
order(global) order of the method (corresponding to b)
apointer to s*s real functions interpreted as a_{ij}=a[i*s+j]
bpointer to s real functions interpreted as b_j=b[j]
btpointer to s real functions interpreted as bt_j=bt[j]
cpointer to s real numbers interpreted as c_i=c[i]

Member Function Documentation

◆ a()

template<class real_type >
function_type dg::mat::FunctionalButcherTableau< real_type >::a ( unsigned  i,
unsigned  j 
) const
inline

Read the a_ij coefficients.

Parameters
irow number 0<=i<s, i>=s results in undefined behaviour
jcol 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 >
function_type dg::mat::FunctionalButcherTableau< real_type >::b ( unsigned  j) const
inline

Read the b_j coefficients.

Parameters
jcol 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 >
function_type dg::mat::FunctionalButcherTableau< real_type >::bt ( unsigned  j) const
inline

Read the embedded bt_j coefficients.

Parameters
jcol 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 >
real_type dg::mat::FunctionalButcherTableau< real_type >::c ( unsigned  i) const
inline

Read the c_i coefficients.

Parameters
irow number 0<=i<s, i>=s results in undefined behaviour
Returns
c_i

◆ embedded_order()

template<class real_type >
unsigned dg::mat::FunctionalButcherTableau< real_type >::embedded_order ( ) const
inline

global order of accuracy for the embedded method represented by bt

◆ isEmbedded()

template<class real_type >
bool dg::mat::FunctionalButcherTableau< real_type >::isEmbedded ( ) const
inline

True if the method has an embedding.

◆ isImplicit()

template<class real_type >
bool dg::mat::FunctionalButcherTableau< real_type >::isImplicit ( ) const
inline

True if an element on or above the diagonal in a is non-zero.

◆ num_stages()

template<class real_type >
unsigned dg::mat::FunctionalButcherTableau< real_type >::num_stages ( ) const
inline

The number of stages s.

◆ order()

template<class real_type >
unsigned dg::mat::FunctionalButcherTableau< real_type >::order ( ) const
inline

global order of accuracy for the method represented by b


The documentation for this struct was generated from the following file: