ampsci
c++ program for high-precision atomic structure calculations of single-valence systems
Public Member Functions | List of all members
AdamsMoulton::DerivativeMatrix< T, Y >abstract

Pure-virtual struct, holds the derivative matrix for 2x2 system of ODEs. Derive from this, and implement a(t),b(t),c(t),d(t) to define the 2x2 ODE. More...

#include <AdamsMoulton.hpp>

Public Member Functions

virtual Y a (T t) const =0
 a,b,c,d are derivative matrix functions; all must be user implemented
 
virtual Y b (T t) const =0
 
virtual Y c (T t) const =0
 
virtual Y d (T t) const =0
 
virtual Y Sf (T) const
 Sf and Sg are optional inhomogenous terms.
 
virtual Y Sg (T) const
 

Detailed Description

template<typename T = double, typename Y = double>
struct AdamsMoulton::DerivativeMatrix< T, Y >

Pure-virtual struct, holds the derivative matrix for 2x2 system of ODEs. Derive from this, and implement a(t),b(t),c(t),d(t) to define the 2x2 ODE.

This struct is used by ODESolver2D, and defined the Derivative matrix, and the (optional) inhomogenous term.

The system of ODEs is defined such that:

\[ \frac{dF(t)}{dt} = D(t) * F(t) + S(t) \]

Where F is a 2D set of functions:

\[ F(t) = \begin{pmatrix} f(t)\\ g(t) \end{pmatrix}, \]

D is the 2x2 "derivative matrix":

\[ D(t) = \begin{pmatrix} a(t) & b(t)\\ c(t) & d(t) \end{pmatrix}, \]

and S(t) is the (optional) 2D inhomogenous term:

\[ S(t) = \begin{pmatrix} s_f(t)\\ s_g(t) \end{pmatrix}. \]

D (and, optionally S) must be provided by the user by implementing this DerivativeMatrix.
D is defined by four functions, a,b,c,d. These four functions must be overriden with definitions to define the ODE system.

Template parameters, T and Y

Each of the four functions is a function of t, which has type T. Usually, T=double or complex<double>, but it may also be an index type (e.g., std::size_t) if the derivative matrix is only known numerically at certain grid points/stored in an array.

The return type, Y, is also a template parameter. It is usually double, but may be any type (e.g., float or complex<double>). DerivativeMatrix will work with any type, though only arithmetic types will work in the ODE solver.

See documentation of ODESolver2D for example.

Note: in tests, deriving from a struct was significantly faster than using std::function, slightly faster than using function pointers, and performed about equally to directly implementing the DerivativeMatrix.


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