|
libode
Easy to compile, fast ODE integrators as C++ classes
|
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 | |
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.
| 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:
| [in] | U | n by n upper triangular matrix |
| [in] | b | right hand side of matrix equation U*out = b |
| [in] | n | dimensions of A and length of b |
| [out] | out | array of length n, result of substitution |
Definition at line 18 of file ode_linalg.cc.
| void ode::ode_check_write | ( | const char * | fn | ) |
| 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.
| [in] | L | n by n lower triangular matrix |
| [in] | b | right hand side of matrix equation U*out = b |
| [in] | p | array of permutation indices from crout_LU |
| [in] | n | dimensions of A and length of b |
| [out] | out | array of length n, result of substitution |
Definition at line 7 of file ode_linalg.cc.
| 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:
| A | matrix of size n to decompose (overwritten) | |
| [in] | n | size of square matrix |
| [out] | p | integer array loaded with permutation indices |
Definition at line 29 of file ode_linalg.cc.
| std::string ode::ode_int_to_string | ( | long | i | ) |
| 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.
| double ode::ode_max2 | ( | double | a, |
| double | b ) |
Simple maximum of two doubles.
Definition at line 7 of file ode_util.cc.
| double ode::ode_min2 | ( | double | a, |
| double | b ) |
Simple minimum of two doubles.
Definition at line 12 of file ode_util.cc.
| void ode::ode_print_exit | ( | const char * | msg | ) |
| void ode::ode_solve_A | ( | double ** | A, |
| double * | b, | ||
| int | n, | ||
| double * | out ) |
Solves a matrix equation A x = b with crout LU decomposition.
| [in] | A | matrix, is overwritten |
| [in] | b | right-hand side of matrix equation |
| [in] | n | size of matrix and arrays |
| [out] | out | solution to linear system |
Definition at line 84 of file ode_linalg.cc.
| 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.
| [in] | LU | crout LU decomposed matrix |
| [in] | p | permutation indices |
| [in] | b | right-hand side of matrix equation |
| [in] | n | size of matrix and arrays |
| [out] | out | solution to linear system |
Definition at line 76 of file ode_linalg.cc.
| 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.
| [in] | T | tridiagonal matrix in banded storage form |
| [in] | r | right-hand side of matrix equation |
| [in] | temp | temporary array |
| [in] | n | size of matrix and arrays |
| [out] | out | solution to linear system |
Definition at line 98 of file ode_linalg.cc.