Libode
High-order ODE solvers as C++ classes with no dependencies
Install / Use
/learn @markmbaum/LibodeREADME
libode
| Reference | Status | Docs |
| --------- | ------ | ---- |
| |
|
|
Easy to compile ODE integrators as C++ classes
This repo contains a collection of C++ classes for solving systems of ordinary differential equations (ODEs) in autonomous form. Documentation is here. All of the solvers are single-step, Runge-Kutta-like methods. There are explicit, adaptive solvers up to the ninth order. The repository also includes Rosenbrock methods, a singly-diagonal implicit Runge-Kutta (SDIRK) method, and several fully implicit Runge-Kutta methods. However, only a few of the implicit methods have solid adaptive time steppers at this point. With the current collection of solvers and features, libode is well suited to any non-stiff systems and to stiff systems that are tightly coupled and have a known Jacobian (ones that don't require sparse or banded matrix routines). It's been useful for solving the same system a huge number of times with varying parameters, when the speed advantage of parallel C++ might be worth it.
The classes were originally styled after Chris Rycroft's example classes. Their structure makes it easy to build a templated integrator on top of an arbitrary solver class and switch the solver/method. Implicit methods can be given a function for the ODE system's Jacobian or, if none is provided, the Jacobian is estimated using finite differences.
Several of the solvers and lots of details on the methods can be found in these amazing books:
- Hairer, E., Nørsett, S. P. & Wanner, G. Solving Ordinary Differential Equations I: Nonstiff Problems. (Springer-Verlag, 1987).
- Hairer, E. & Wanner, G. Solving Ordinary Differential Equations II: Stiff and Differential-Algebraic Problems. (Springer, 1996).
The table below lists all the solvers and gives some basic information about them. All of the solvers can be used with a customized time step selection function, but those with a built-in adaptive capability are indicated below. Papers and/or links to the derivation or original publication of the methods are often copied in the headers for the solver classes and included in the documentation. Some work still needs to be done to make the implicit methods genuinely useful, and a list of things to implement is in the todo.md file.
Method | Class Name | Header File | (ex/im)plicit | built-in adaptive? | stages | order | stability
--- | --- | --- | --- | --- | --- | --- | ---
Forward Euler | OdeEuler | ode_euler.h | explicit | no | 1 | 1
Trapezoidal Rule | OdeTrapz | ode_trapz.h | explicit | no | 2 | 2
Strong Stability-Preserving, Order 3 | OdeSsp3 | ode_ssp_3.h | explicit | no | 3 | 3
Runge-Kutta-Fehlberg 3(2) | OdeRKF32| ode_rkf_32.h | explicit | yes | 3 | 3
RK4 | OdeRK4 | ode_rk_4.h | explicit | no | 4 | 4
Runge-Kutta 4(3) | OdeRK43| ode_rk_43.h | explicit | yes | 5 | 4
Cash-Karp | OdeRKCK | ode_rkck.h | explicit | yes | 6 | 5
Dormand-Prince 5(4) | OdeDoPri54| ode_dopri_54.h | explicit | yes | 7 | 5
Jim Verner's "most efficent" 6(5) | OdeVern65 | ode_vern_65.h | explicit | yes | 9 | 6
Jim Verner's "most efficent" 7(6) | OdeVern76 | ode_vern_76.h | explicit | yes | 10 | 7
Dormand-Prince 8(7) | OdeDoPri87 | ode_dopri_87.h | explicit | yes | 13 | 8
Jim Verner's "most efficent" 9(8) | OdeVern98 | ode_vern_98.h | explicit | yes | 16 | 9
Rosenbrock 4(3) | OdeGRK4A | ode_grk4a.h | implicit | yes | 4 | 4 | A
Rosenbrock 6 | OdeROW6A | ode_row6a.h | implicit | no | 6 | 6 | A
Backward Euler | OdeBackwardEuler | ode_backward_euler.h | implicit | no | 1 | 1 | L
Gauss 6th Order | OdeGauss6 | ode_gauss_6.h | implicit | not yet | 3 | 6 | A
Lobatto IIIC 6th Order | OdeLobattoIIIC6 | ode_lobatto_iiic_6.h | implicit | not yet | 4 | 6 | L
Radau IIA 5th Order | OdeRadauIIA5 | ode_radau_iia_5.h | implicit | not yet | 3 | 5 | L
Geng's Symplectic 5th Order | OdeGeng5 | ode_geng_5.h | implicit | no | 3 | 5 | A?
SDIRK 4(3) | OdeSDIRK43 | ode_sdirk_43.h | implicit | yes | 4 | 4 | L
Compiling
Compiling with Makefile
Short Instructions
- Set your compiler and flags by defining
CXXandCFLAGSenvironment variables or by uncommenting and editing those variables in the config.mk file. - Run
makein the top directory where the Makefile is. - If anything weird happens, tell me.
- Run
make testsand execute therun_tests.shscript to check that things are working. - If you want, also execute
run_examples.shto run some example solvers (Python with numpy and matplotlib are needed for plotting). To do this, you must have theCXXandCFLAGSvariables active in the config.mk file. - Create derived classes and link to the library with
-I<path>/libode/include/ode -L<path>/libode/bin -lode, replacing<path>with the path to the directory abovelibodeon your computer.
Longer Instructions
libode is meant to provide straightforward access to class-based ODE solvers without dependencies or specialized compiling processes. The library is free-standing and there is only one step to take before simply running the Makefile and being done with it. Consequently, the library is also slim on features and doesn't provide things like sparse matrices and dense output. For many systems of ODEs, though, libode should make it easy to build an integrator and enjoy the speed of C++ and openmp without the headaches of large, complex packages.
First, before any of the libode classes can be compiled, you must tell the Makefile which compiler and compiler flags to use. This can be done two ways. First, you can set the environment variables CXX and CFLAGS. For example, in the shell:
CXX=g++
CFLAGS=-Wall -O3
Then the Makefile will use the g++ compiler and always include the -Wall and -O3 flags. Second, you can uncomment and edit variables in the config.mk file. Editing config.mk is required if you want to run the run_examples.sh script, which solves and plots some example systems of ODEs.
Then simply run make and everything should compile.
The Makefile compiles all of the necessary code into the obj folder, then archives it in the bin directory as a file called libode.a. To use the solvers, you can link libode.a or the object files directly (in the obj directory) when compiling your derived class. You must also include the appropriate header files from the src directory, as there is no unified header file for the library. All of the classes have their header file name displayed in the documentation and in the table above. Linking the solver classes requires something like
-I<path>/libode/include/ode -L<path>/libode/bin -lode
when compiling derived code, with <path> replaced by path elements leading to the libode directory. For some examples of how to link a derived class to libode and create a program to run integrations, see the examples folder.
Test programs are compiled with make tests and they can all be run in sequence with the run_tests.sh script.
Compiling with CMake
Whe library can be built with CMake as well. You can start CMake with the following script
. run_cmake.sh
This script creates an folder build, compiles the library, installs it locally and creates a package.
Here is an example of an file CMakeLists.txt to show how to use the ode library:
cmake_minimum_required(VERSION 3.15)
project("ODETest" VERSION 0.9 DESCRIPTION "A project with the libode library")
set(CMAKE_CXX_STANDARD 11)
set(CMAKE_CXX_STANDARD_REQUIRED ON)
# It is a local prefix.
# Comment it out if you want to use the default prefix.
set(CMAKE_PREFIX_PATH "$ENV{HOME}/install")
find_package(ode CONFIG REQUIRED)
add_executable(${PROJECT_NAME})
target_sources(${PROJECT_NAME}
PRIVATE
my_ode_test.cpp # These are the names of the source files.
)
target_link_libraries(${PROJECT_NAME} PRIVATE ode::ode)
Examples
Several example programs for interesting/famous systems of ODEs are in the examples folder. In each of the example directories, the programs can be compiled, executed, and plotted simply by running the run.sh script (assuming the config.mk file is set up and an up-to-date version of Python is installed on your computer with numpy, scipy, and matplotlib). These programs are good examples of how to put everything together and use the solvers. To run all the examples in sequence and look at the plotted results, run the run_examples.sh script.
Using the Solvers
Define a Class
To integrate a specific system of ODEs, a new class must be created to inherit from one of the solver classes. This new inheriting class must
- Define the system of ODEs to be solved by implementing the
ode_fun()function. This is a virtual function in the base classes. Once implemented, it's used by the stepping and solving functions. - Set initial conditions using the
set_sol()function. This can be done inside the new class's constructor, or not. - Optionally implement the
ode_jac()function for implicit methods. This is also a virtual function in the base classes. If it's not overridden but is needed, a (crude) finite-difference estimate of the Jacobian is used.
For flexibility, the derived class can be a template, so that the solver/method can be chosen when the class is constructed. Other than defining the system of equations and setting initial conditions, the derived class can store whatever informatio
