Discontinuous Galerkin Library
#include "dg/algorithm.h"
|
Runge-Kutta explicit ODE-integrators. More...
#include <cassert>
#include <array>
#include <tuple>
#include "ode.h"
#include "backend/exceptions.h"
#include "tableau.h"
#include "blas1.h"
#include "implicit.h"
Go to the source code of this file.
Classes | |
struct | dg::IdentityFilter |
A filter that does nothing. More... | |
struct | dg::ERKStep< ContainerType > |
Embedded Runge Kutta explicit time-step with error estimate \(
\begin{align}
k_i = f\left( t^n + c_i \Delta t, u^n + \Delta t \sum_{j=1}^{i-1} a_{ij} k_j\right) \\
u^{n+1} = u^{n} + \Delta t\sum_{j=1}^s b_j k_j \\
\tilde u^{n+1} = u^{n} + \Delta t\sum_{j=1}^s \tilde b_j k_j \\
\delta^{n+1} = u^{n+1} - \tilde u^{n+1} = \Delta t\sum_{j=1}^s (b_j - \tilde b_j) k_j
\end{align}
\). More... | |
struct | dg::FilteredERKStep< ContainerType > |
EXPERIMENTAL: Filtered Embedded Runge Kutta explicit time-step with error estimate \(
\begin{align}
k_i = f\left( t^n + c_i \Delta t, \Lambda\Pi \left[u^n + \Delta t \sum_{j=1}^{i-1} a_{ij} k_j\right]\right) \\
u^{n+1} = \Lambda\Pi\left[u^{n} + \Delta t\sum_{j=1}^s b_j k_j\right] \\
\delta^{n+1} = \Delta t\sum_{j=1}^s (\tilde b_j - b_j) k_j
\end{align}
\). More... | |
struct | dg::ARKStep< ContainerType > |
Additive Runge Kutta (semi-implicit) time-step with error estimate following The ARKode library More... | |
struct | dg::DIRKStep< ContainerType > |
Embedded diagonally implicit Runge Kutta time-step with error estimate \(
\begin{align}
k_i = f\left( t^n + c_i \Delta t, u^n + \Delta t \sum_{j=1}^{i} a_{ij} k_j\right) \\
u^{n+1} = u^{n} + \Delta t\sum_{j=1}^s b_j k_j \\
\tilde u^{n+1} = u^{n} + \Delta t\sum_{j=1}^s \tilde b_j k_j \\
\delta^{n+1} = u^{n+1} - \tilde u^{n+1} = \Delta t\sum_{j=1}^s (b_j-\tilde b_j) k_j
\end{align}
\). More... | |
struct | dg::ShuOsher< ContainerType > |
Shu-Osher fixed-step explicit ODE integrator with Slope Limiter / Filter \(
\begin{align}
u_0 &= u_n \\
u_i &= \Lambda\Pi \left(\sum_{j=0}^{i-1}\left[ \alpha_{ij} u_j + \Delta t \beta_{ij} f( t_j, u_j)\right]\right)\\
u^{n+1} &= u_s
\end{align}
\). More... | |
struct | dg::SinglestepTimeloop< ContainerType > |
Integrate using a for loop and a fixed time-step. More... | |
Namespaces | |
namespace | dg |
This is the namespace for all functions and classes defined and used by the discontinuous Galerkin library. | |
Typedefs | |
template<class ContainerType > | |
using | dg::RungeKutta = ERKStep<ContainerType> |
Runge-Kutta fixed-step explicit ODE integrator \(
\begin{align}
k_i = f\left( t^n + c_i \Delta t, u^n + \Delta t \sum_{j=1}^{i-1} a_{ij} k_j\right) \\
u^{n+1} = u^{n} + \Delta t\sum_{j=1}^s b_j k_j
\end{align}
\). | |
template<class ContainerType > | |
using | dg::FilteredRungeKutta = FilteredERKStep<ContainerType> |
Filtered Runge-Kutta fixed-step explicit ODE integrator \(
\begin{align}
k_i = f\left( t^n + c_i \Delta t, \Lambda\Pi \left[u^n + \Delta t \sum_{j=1}^{i-1} a_{ij} k_j\right]\right) \\
u^{n+1} = \Lambda\Pi\left[u^{n} + \Delta t\sum_{j=1}^s b_j k_j\right]
\end{align}
\). | |
template<class ContainerType > | |
using | dg::ImplicitRungeKutta = DIRKStep<ContainerType> |
Runge-Kutta fixed-step implicit ODE integrator \(
\begin{align}
k_i = f\left( t^n + c_i \Delta t, u^n + \Delta t \sum_{j=1}^{s} a_{ij} k_j\right) \\
u^{n+1} = u^{n} + \Delta t\sum_{j=1}^s b_j k_j
\end{align}
\). | |
Functions | |
bool | dg::is_same (double x, double y, double eps=1e-15) |
bool | dg::is_same (float x, float y, float eps=1e-6) |
template<class T > | |
bool | dg::is_same (T x, T y) |
bool | dg::is_divisable (double a, double b, double eps=1e-15) |
bool | dg::is_divisable (float a, float b, float eps=1e-6) |
Runge-Kutta explicit ODE-integrators.