0.2.0
Loading...
Searching...
No Matches
NumericalIntegration.hpp File Reference

Provides Numerical integration methods. More...

Go to the source code of this file.

Functions

template<typename Y , typename Z , typename Scalar , typename = std::enable_if_t<std::is_floating_point_v<Scalar>>>
NAV::Heun2 (const Y &y_n, const std::array< Z, 2 > z, const Scalar &h, const auto &f, const auto &constParam, const Scalar &t_n=0)
 Heun's method (2nd order) (explicit) .
 
template<typename Y , typename Z , typename Scalar , typename = std::enable_if_t<std::is_floating_point_v<Scalar>>>
NAV::Heun3 (const Y &y_n, const std::array< Z, 3 > z, const Scalar &h, const auto &f, const auto &constParam, const Scalar &t_n=0)
 Heun's method (3nd order) (explicit) .
 
template<typename Y , typename Z , typename Scalar , typename = std::enable_if_t<std::is_floating_point_v<Scalar>>>
NAV::RungeKutta1 (const Y &y_n, const std::array< Z, 1 > z, const Scalar &h, const auto &f, const auto &constParam, const Scalar &t_n=0)
 Runge-Kutta 1st order (explicit) / (Forward) Euler method .
 
template<typename Y , typename Z , typename Scalar , typename = std::enable_if_t<std::is_floating_point_v<Scalar>>>
NAV::RungeKutta2 (const Y &y_n, const std::array< Z, 2 > z, const Scalar &h, const auto &f, const auto &constParam, const Scalar &t_n=0)
 Runge-Kutta 2nd order (explicit) / Explicit midpoint method .
 
template<typename Y , typename Z , typename Scalar , typename = std::enable_if_t<std::is_floating_point_v<Scalar>>>
NAV::RungeKutta3 (const Y &y_n, const std::array< Z, 3 > z, const Scalar &h, const auto &f, const auto &constParam, const Scalar &t_n=0)
 Runge-Kutta 3rd order (explicit) / Simpson's rule .
 
template<typename Y , typename Z , typename Scalar , typename = std::enable_if_t<std::is_floating_point_v<Scalar>>>
NAV::RungeKutta4 (const Y &y_n, const std::array< Z, 4 > z, const Scalar &h, const auto &f, const auto &constParam, const Scalar &t_n=0)
 Runge-Kutta 4th order (explicit) .
 
template<typename Y , typename Z , typename Scalar , size_t s, std::array< std::array< Scalar, s+1 >, s+1 > butcherTableau, typename = std::enable_if_t<std::is_floating_point_v<Scalar>>>
NAV::ButcherTableau::RungeKuttaExplicit (const Y &y_n, const std::array< Z, s > z, const Scalar &h, const auto &f, const auto &constParam, const Scalar &t_n)
 Calculates explicit Runge-Kutta methods. Order is defined by the Butcher tableau.
 

Variables

template<typename Scalar >
constexpr std::array< std::array< Scalar, 3 >, 3 > NAV::ButcherTableau::Heun2
 Butcher tableau for Heun's method (2nd order) (explicit)
 
template<typename Scalar >
constexpr std::array< std::array< Scalar, 4 >, 4 > NAV::ButcherTableau::Heun3
 Butcher tableau for Heun's method (3nd order) (explicit)
 
template<typename Scalar >
constexpr std::array< std::array< Scalar, 2 >, 2 > NAV::ButcherTableau::RK1
 Butcher tableau for Runge-Kutta 1st order (explicit) / (Forward) Euler method.
 
template<typename Scalar >
constexpr std::array< std::array< Scalar, 3 >, 3 > NAV::ButcherTableau::RK2
 Butcher tableau for Runge-Kutta 2nd order (explicit) / Explicit midpoint method.
 
template<typename Scalar >
constexpr std::array< std::array< Scalar, 4 >, 4 > NAV::ButcherTableau::RK3
 Butcher tableau for Runge-Kutta 3rd order (explicit) / Simpson's rule.
 
template<typename Scalar >
constexpr std::array< std::array< Scalar, 5 >, 5 > NAV::ButcherTableau::RK4
 Butcher tableau for Runge-Kutta 4th order (explicit)
 

Detailed Description

Provides Numerical integration methods.

Author
T. Topp (topp@.nosp@m.ins..nosp@m.uni-s.nosp@m.tutt.nosp@m.gart..nosp@m.de)
Date
2023-12-11

Function Documentation

◆ Heun2()

template<typename Y , typename Z , typename Scalar , typename = std::enable_if_t<std::is_floating_point_v<Scalar>>>
Y NAV::Heun2 ( const Y & y_n,
const std::array< Z, 2 > z,
const Scalar & h,
const auto & f,
const auto & constParam,
const Scalar & t_n = 0 )

Heun's method (2nd order) (explicit) .

\begin{equation} \label{eq:eq-Heun2} \begin{aligned} y_{n+1} &= y_n + \frac{h}{2} (k_1 + k_2) \\[1em] k_1 &= f(t_n, y_n) \\ k_2 &= f(t_n + h, y_n + h k_1) \\ \end{aligned} \end{equation}

Butcher tableau:

\begin{equation} \label{eq:eq-Heun2-bt} \renewcommand\arraystretch{1.2} \begin{array}{c|cc} 0 \\ 1 & 1 \\ \hline & \frac{1}{2} & \frac{1}{2} \\ \end{array} \end{equation}

Parameters
[in]y_nState vector at time t_n
[in]zArray of measurements, one for each evaluation point of the Runge Kutta
[in]hIntegration step in [s]
[in]fTime derivative function
[in]constParamConstant parameters passed to each time derivative function call
[in]t_nTime t_n

◆ Heun3()

template<typename Y , typename Z , typename Scalar , typename = std::enable_if_t<std::is_floating_point_v<Scalar>>>
Y NAV::Heun3 ( const Y & y_n,
const std::array< Z, 3 > z,
const Scalar & h,
const auto & f,
const auto & constParam,
const Scalar & t_n = 0 )

Heun's method (3nd order) (explicit) .

\begin{equation} \label{eq:eq-Heun3} \begin{aligned} y_{n+1} &= y_n + \frac{h}{4} (k_1 + 3 k_3) \\[1em] k_1 &= f(t_n, y_n) \\ k_2 &= f(t_n + \frac{h}{3}, y_n + \frac{h}{3} k_1) \\ k_3 &= f(t_n + \frac{2 h}{3}, y_n + \frac{2 h}{3} k_2) \\ \end{aligned} \end{equation}

Butcher tableau:

\begin{equation} \label{eq:eq-Heun3-bt} \renewcommand\arraystretch{1.2} \begin{array}{c|ccc} 0 \\ \frac{1}{3} & \frac{1}{3} \\ \frac{2}{3} & 0 & \frac{2}{3} \\ \hline & \frac{1}{4} & 0 & \frac{3}{4} \\ \end{array} \end{equation}

Parameters
[in]y_nState vector at time t_n
[in]zArray of measurements, one for each evaluation point of the Runge Kutta
[in]hIntegration step in [s]
[in]fTime derivative function
[in]constParamConstant parameters passed to each time derivative function call
[in]t_nTime t_n

◆ RungeKutta1()

template<typename Y , typename Z , typename Scalar , typename = std::enable_if_t<std::is_floating_point_v<Scalar>>>
Y NAV::RungeKutta1 ( const Y & y_n,
const std::array< Z, 1 > z,
const Scalar & h,
const auto & f,
const auto & constParam,
const Scalar & t_n = 0 )

Runge-Kutta 1st order (explicit) / (Forward) Euler method .

\begin{equation} \label{eq:eq-RungeKutta1-explicit} y_{n+1} = y_n + h f(t_n, y_n) \end{equation}

Butcher tableau:

\begin{equation} \label{eq:eq-RungeKutta1-explicit-bt} \renewcommand\arraystretch{1.2} \begin{array}{c|c} 0 \\ \hline & 1 \\ \end{array} \end{equation}

Parameters
[in]y_nState vector at time t_n
[in]zArray of measurements, one for each evaluation point of the Runge Kutta
[in]hIntegration step in [s]
[in]fTime derivative function
[in]constParamConstant parameters passed to each time derivative function call
[in]t_nTime t_n

◆ RungeKutta2()

template<typename Y , typename Z , typename Scalar , typename = std::enable_if_t<std::is_floating_point_v<Scalar>>>
Y NAV::RungeKutta2 ( const Y & y_n,
const std::array< Z, 2 > z,
const Scalar & h,
const auto & f,
const auto & constParam,
const Scalar & t_n = 0 )

Runge-Kutta 2nd order (explicit) / Explicit midpoint method .

\begin{equation} \label{eq:eq-RungeKutta2-explicit} \begin{aligned} y_{n+1} &= y_n + h k_2 \\[1em] k_1 &= f(t_n, y_n) \\ k_2 &= f(t_n + \frac{h}{2}, y_n + \frac{h}{2} k_1) \\ \end{aligned} \end{equation}

Butcher tableau:

\begin{equation} \label{eq:eq-RungeKutta2-explicit-bt} \renewcommand\arraystretch{1.2} \begin{array}{c|cc} 0 \\ \frac{1}{2} & \frac{1}{2} \\ \hline & 0 & 1 \\ \end{array} \end{equation}

Parameters
[in]y_nState vector at time t_n
[in]zArray of measurements, one for each evaluation point of the Runge Kutta
[in]hIntegration step in [s]
[in]fTime derivative function
[in]constParamConstant parameters passed to each time derivative function call
[in]t_nTime t_n

◆ RungeKutta3()

template<typename Y , typename Z , typename Scalar , typename = std::enable_if_t<std::is_floating_point_v<Scalar>>>
Y NAV::RungeKutta3 ( const Y & y_n,
const std::array< Z, 3 > z,
const Scalar & h,
const auto & f,
const auto & constParam,
const Scalar & t_n = 0 )

Runge-Kutta 3rd order (explicit) / Simpson's rule .

\begin{equation} \label{eq:eq-RungeKutta3-explicit} \begin{aligned} y_{n+1} &= y_n + \frac{h}{6} ( k_1 + 4 k_2 + k_3 ) \\[1em] k_1 &= f(t_n, y_n) \\ k_2 &= f(t_n + \frac{h}{2}, y_n + \frac{h}{2} k_1) \\ k_3 &= f(t_n + h, y_n + h (-k_1 + 2 k_2)) \\ \end{aligned} \end{equation}

Butcher tableau:

\begin{equation} \label{eq:eq-RungeKutta3-explicit-bt} \renewcommand\arraystretch{1.2} \begin{array}{c|ccc} 0 \\ \frac{1}{2} & \frac{1}{2} \\ 1 & -1 & 2 \\ \hline & \frac{1}{6} & \frac{4}{6} & \frac{1}{6} \\ \end{array} \end{equation}

Parameters
[in]y_nState vector at time t_n
[in]zArray of measurements, one for each evaluation point of the Runge Kutta
[in]hIntegration step in [s]
[in]fTime derivative function
[in]constParamConstant parameters passed to each time derivative function call
[in]t_nTime t_n

◆ RungeKutta4()

template<typename Y , typename Z , typename Scalar , typename = std::enable_if_t<std::is_floating_point_v<Scalar>>>
Y NAV::RungeKutta4 ( const Y & y_n,
const std::array< Z, 4 > z,
const Scalar & h,
const auto & f,
const auto & constParam,
const Scalar & t_n = 0 )

Runge-Kutta 4th order (explicit) .

\begin{equation} \label{eq:eq-RungeKutta4-explicit} \begin{aligned} y_{n+1} &= y_n + \frac{h}{6} ( k_1 + 2 k_2 + 2 k_3 + k_4 ) \\[1em] k_1 &= f(t_n, y_n) \\ k_2 &= f(t_n + \frac{h}{2}, y_n + \frac{h}{2} k_1) \\ k_3 &= f(t_n + \frac{h}{2}, y_n + \frac{h}{2} k_2) \\ k_4 &= f(t_n + h, y_n + h k_3) \\ \end{aligned} \end{equation}

Butcher tableau:

\begin{equation} \label{eq:eq-RungeKutta4-explicit-bt} \renewcommand\arraystretch{1.2} \begin{array}{c|cccc} 0 \\ \frac{1}{2} & \frac{1}{2} \\ \frac{1}{2} & 0 & \frac{1}{2} \\ 1 & 0 & 0 & 1 \\ \hline & \frac{1}{6} & \frac{1}{3} & \frac{1}{3} & \frac{1}{6} \\ \end{array} \end{equation}

Parameters
[in]y_nState vector at time t_n
[in]zArray of measurements, one for each evaluation point of the Runge Kutta
[in]hIntegration step in [s]
[in]fTime derivative function
[in]constParamConstant parameters passed to each time derivative function call
[in]t_nTime t_n

◆ RungeKuttaExplicit()

template<typename Y , typename Z , typename Scalar , size_t s, std::array< std::array< Scalar, s+1 >, s+1 > butcherTableau, typename = std::enable_if_t<std::is_floating_point_v<Scalar>>>
Y NAV::ButcherTableau::RungeKuttaExplicit ( const Y & y_n,
const std::array< Z, s > z,
const Scalar & h,
const auto & f,
const auto & constParam,
const Scalar & t_n )
inline

Calculates explicit Runge-Kutta methods. Order is defined by the Butcher tableau.

Parameters
[in]y_nState vector at time t_n
[in]zArray of measurements, one for each evaluation point of the Runge Kutta
[in]hIntegration step in [s]
[in]fTime derivative function
[in]constParamConstant parameters passed to each time derivative function call
[in]t_nTime t_n
Returns
State vector at time t_(n+1)