libode
Easy to compile, fast ODE integrators as C++ classes
Loading...
Searching...
No Matches
ode Namespace Reference

Classes

class  NewtonBackwardEuler
 Nonlinear system solver for OdeBackwardEuler. More...
 
class  NewtonGauss6
 Nonlinear system solver for OdeGauss6. More...
 
class  NewtonGeng5
 Nonlinear system solver for OdeGeng5. More...
 
class  NewtonLobattoIIIC6
 Nonlinear system solver for OdeLobattoIIIC6. More...
 
class  NewtonRadauIIA5
 Nonlinear system solver for OdeRadauIIA5. More...
 
class  NewtonSDIRK43
 Nonlinear system solver for OdeSDIRK43. More...
 
class  OdeAdaptive
 Base class implementing solver functions with adaptive time steps. More...
 
class  OdeBackwardEuler
 Backward Euler's method, unconditionally stable but relatively inaccurate. More...
 
class  OdeBase
 Lowest base class for all solvers. More...
 
class  OdeDoPri54
 Popular explicit 5/4 pair from Dormand & Prince. More...
 
class  OdeDoPri87
 Explicit 8/7 pair from Dormand & Prince. More...
 
class  OdeEmbedded
 Base clase implementing methods for embedded Runge-Kutta error estimation. More...
 
class  OdeERK
 Base class providing space for temporary solutions moving through RK stages. More...
 
class  OdeEuler
 The simplest runge kutta method, forward Euler's. More...
 
class  OdeGauss6
 The sixth-order, A-stable, fully-implicit Gauss-Legendre method with 3 stages. More...
 
class  OdeGeng5
 The fifth-order, symplectic, fully-implicit Geng integrator with 3 stages. More...
 
class  OdeGRK4A
 Fourth-order, A-stable, adaptive Rosenbrock method from Kaps and Rentrop. More...
 
class  OdeIRK
 Provides a large vector containing the slope values of all stages with pointers to each of the individual stages. More...
 
class  OdeLobattoIIIC6
 The sixth-order, L-stable, fully-implicit Lobatto IIIC method with 4 stages. More...
 
class  OdeNewton
 Newton's method for nonlinear systems of equations. More...
 
class  OdeNewtonBridge
 Templated base class connecting solver objects and OdeNewton objects. More...
 
class  OdeNewtonIRK
 Extension of OdeNewtonBridge class for fully implicit methods. More...
 
class  OdeNewtonSDIRK
 Extension of OdeNewtonBridge class for fully SDIRK methods. More...
 
class  OdeRadauIIA5
 The fifth-order, L-stable, fully-implicit Radau IIA method with 3 stages. More...
 
class  OdeRK
 Provides space for stage slope values, an array of arrays for k values. More...
 
class  OdeRK4
 The classic Runge-Kutta 4th order method. More...
 
class  OdeRK43
 This class implements a 3rd and 4th order method with the FSAL (first same as last) property. More...
 
class  OdeRKCK
 Explicit 5/4 pair, also with 3rd, 2nd, and 1st order embedded methods, from Cash & Karp. More...
 
class  OdeRKF32
 2nd and 3rd order solver developed by Fehlberg More...
 
class  OdeRosenbrock
 Base class for Rosenbrock methods. More...
 
class  OdeROW6A
 6th order, A-stable Rosenbrock method More...
 
class  OdeSDIRK43
 L-stable 4/3 SDIRK pair. More...
 
class  OdeSsp3
 Strong stability preserving method of order 3. More...
 
class  OdeTrapz
 Second order, explicit trapezoidal rule. More...
 
class  OdeVern65
 Jim Verner's "most efficient" 6/5 pair. More...
 
class  OdeVern76
 Jim Verner's "most efficient" 7/6 pair. More...
 
class  OdeVern98
 Jim Verner's "most efficient" 9/8 pair. More...
 

Functions

void ode_check_write (const char *fn)
 check if a file can be written (opened in write mode)
 
void ode_print_exit (const char *msg)
 print a message and exit with failure
 
std::string ode_int_to_string (long i)
 converts an integer into a string
 
void ode_crout_forw_sub (double **L, double *b, int *p, int n, double *out)
 Crout forward substitution.
 
void ode_back_sub (double **U, double *b, int n, double *out)
 Backward substitution.
 
void ode_crout_LU (double **A, int n, int *p)
 Crout LU decomposition.
 
void ode_solve_LU (double **LU, int *p, double *b, int n, double *out)
 Solves a matrix equation where the matrix has already be crout LU decomposed.
 
void ode_solve_A (double **A, double *b, int n, double *out)
 Solves a matrix equation A x = b with crout LU decomposition.
 
void ode_solve_tridiag (double **T, double *r, double *temp, int n, double *out)
 Solves a tridiagonal matrix equation.
 
double ode_max2 (double a, double b)
 Simple maximum of two doubles.
 
double ode_min2 (double a, double b)
 Simple minimum of two doubles.
 
bool ode_is_close (double a, double b, double thresh)
 Checks if two numbers are very close to each other.
 
template<class T >
void ode_write (char const *fn, T *a, unsigned long n)
 write an array to a binary file
 

Detailed Description

These classes provides some variables of convenience when connecting the Newton base (which is left as a generic Newton solver for potential reuse elsewhere) and the systems required for implicit methods.

Function Documentation

◆ ode_back_sub()

void ode::ode_back_sub ( double ** U,
double * b,
int n,
double * out )

Backward substitution.

Performs backward substitution to solve an upper triangular matrix equation input:

Parameters
[in]Un by n upper triangular matrix
[in]bright hand side of matrix equation U*out = b
[in]ndimensions of A and length of b
[out]outarray of length n, result of substitution

Definition at line 18 of file ode_linalg.cc.

◆ ode_check_write()

void ode::ode_check_write ( const char * fn)

check if a file can be written (opened in write mode)

Parameters
[in]fndesired file name

Definition at line 7 of file ode_io.cc.

◆ ode_crout_forw_sub()

void ode::ode_crout_forw_sub ( double ** L,
double * b,
int * p,
int n,
double * out )

Crout forward substitution.

Performs forward substitution to solve a lower triangular matrix equation where the elements on the main diagonal of the lower triangular matrix are assumed to all be 1.

Parameters
[in]Ln by n lower triangular matrix
[in]bright hand side of matrix equation U*out = b
[in]parray of permutation indices from crout_LU
[in]ndimensions of A and length of b
[out]outarray of length n, result of substitution

Definition at line 7 of file ode_linalg.cc.

◆ ode_crout_LU()

void ode::ode_crout_LU ( double ** A,
int n,
int * p )

Crout LU decomposition.

performs LU decomposition with partial pivoting on the matrix A and sets permutation indices in the array p, following the outline in the book:

  • Demmel, James W. Applied numerical linear algebra. Vol. 56. Siam, 1997.
    Parameters
    Amatrix of size n to decompose (overwritten)
    [in]nsize of square matrix
    [out]pinteger array loaded with permutation indices

Definition at line 29 of file ode_linalg.cc.

◆ ode_int_to_string()

std::string ode::ode_int_to_string ( long i)

converts an integer into a string

Definition at line 25 of file ode_io.cc.

◆ ode_is_close()

bool ode::ode_is_close ( double a,
double b,
double thresh )

Checks if two numbers are very close to each other.

Definition at line 17 of file ode_util.cc.

◆ ode_max2()

double ode::ode_max2 ( double a,
double b )

Simple maximum of two doubles.

Definition at line 7 of file ode_util.cc.

◆ ode_min2()

double ode::ode_min2 ( double a,
double b )

Simple minimum of two doubles.

Definition at line 12 of file ode_util.cc.

◆ ode_print_exit()

void ode::ode_print_exit ( const char * msg)

print a message and exit with failure

Parameters
[in]msgthe message to print before exiting the program

Definition at line 17 of file ode_io.cc.

◆ ode_solve_A()

void ode::ode_solve_A ( double ** A,
double * b,
int n,
double * out )

Solves a matrix equation A x = b with crout LU decomposition.

Parameters
[in]Amatrix, is overwritten
[in]bright-hand side of matrix equation
[in]nsize of matrix and arrays
[out]outsolution to linear system

Definition at line 84 of file ode_linalg.cc.

◆ ode_solve_LU()

void ode::ode_solve_LU ( double ** LU,
int * p,
double * b,
int n,
double * out )

Solves a matrix equation where the matrix has already be crout LU decomposed.

Parameters
[in]LUcrout LU decomposed matrix
[in]ppermutation indices
[in]bright-hand side of matrix equation
[in]nsize of matrix and arrays
[out]outsolution to linear system

Definition at line 76 of file ode_linalg.cc.

◆ ode_solve_tridiag()

void ode::ode_solve_tridiag ( double ** T,
double * r,
double * temp,
int n,
double * out )

Solves a tridiagonal matrix equation.

The input matrix T must be in banded storage. This means that the diagonals are placed on rows so that the trigiagonal matrix A

A = | a00  a01   0    0  |
    | a10  a11  a12   0  |
    |  0   a21  a22  a23 |
    |  0    0   a32  a33 |

is stored in banded form as the array T

T = |  *   a01  a12  a23 |
    | a00  a11  a22  a33 |
    | a10  a21  a31   *  |

and the value of the asterisk elements doesn't matter.

Parameters
[in]Ttridiagonal matrix in banded storage form
[in]rright-hand side of matrix equation
[in]temptemporary array
[in]nsize of matrix and arrays
[out]outsolution to linear system

Definition at line 98 of file ode_linalg.cc.

◆ ode_write()

template<class T >
void ode::ode_write ( char const * fn,
T * a,
unsigned long n )

write an array to a binary file

Parameters
[in]fnfile name or path
[in]aarray of values to write
[in]nsize of array or number of elements to write

Definition at line 27 of file ode_io.h.