...one of the most highly
regarded and expertly designed C++ library projects in the
world.
— Herb Sutter and Andrei
Alexandrescu, C++
Coding Standards
Home | Libraries | People | FAQ | More |
First of all, you have to specify the data type that represents a state
x of your system. Mathematically, this usually is
an n-dimensional vector with real numbers or complex numbers as scalar
objects. For odeint the most natural way is to use vector< double >
or vector< complex< double > >
to represent the system state. However, odeint can deal with other container
types as well, e.g. boost::array< double , N >
, as long as it fulfills some requirements
defined below.
To integrate a differential equation numerically, one also has to define the rhs of the equation x' = f(x). In odeint you supply this function in terms of an object that implements the ()-operator with a certain parameter structure. Hence, the straightforward way would be to just define a function, e.g:
/* The type of container used to hold the state vector */ typedef std::vector< double > state_type; const double gam = 0.15; /* The rhs of x' = f(x) */ void harmonic_oscillator( const state_type &x , state_type &dxdt , const double /* t */ ) { dxdt[0] = x[1]; dxdt[1] = -x[0] - gam*x[1]; }
The parameters of the function must follow the example above where x
is the current state, here a two-component
vector containing position q and momentum p
of the oscillator, dxdt
is the derivative x' and should be filled by the function
with f(x), and t
is the current time. Note that in this example t is
not required to calculate f, however odeint expects
the function signature to have exactly three parameters (there are exception,
discussed later).
A more sophisticated approach is to implement the system as a class where the rhs function is defined as the ()-operator of the class with the same parameter structure as above:
/* The rhs of x' = f(x) defined as a class */ class harm_osc { double m_gam; public: harm_osc( double gam ) : m_gam(gam) { } void operator() ( const state_type &x , state_type &dxdt , const double /* t */ ) { dxdt[0] = x[1]; dxdt[1] = -x[0] - m_gam*x[1]; } };
odeint can deal with instances of such classes instead of pure functions which allows for cleaner code.
Numerical integration works iteratively, that means you start at a state x(t) and perform a time-step of length dt to obtain the approximate state x(t+dt). There exist many different methods to perform such a time-step each of which has a certain order q. If the order of a method is q than it is accurate up to term ~dtq that means the error in x made by such a step is ~dtq+1. odeint provides several steppers of different orders, see Stepper overview.
Some of steppers in the table above are special: Some need the Jacobian of the ODE, others are constructed for special ODE-systems like Hamiltonian systems. We will show typical examples and use-cases in this tutorial and which kind of steppers should be applied.
The basic stepper just performs one time-step and doesn't give you any
information about the error that was made (except that you know it is of
order q+1). Such steppers are used with constant step
size that should be chosen small enough to have reasonable small errors.
However, you should apply some sort of validity check of your results (like
observing conserved quantities) because you have no other control of the
error. The following example defines a basic stepper based on the classical
Runge-Kutta scheme of 4th order. The declaration of the stepper requires
the state type as template parameter. The integration can now be done by
using the integrate_const( Stepper, System, state, start_time, end_time, step_size
)
function from odeint:
runge_kutta4< state_type > stepper; integrate_const( stepper , harmonic_oscillator , x , 0.0 , 10.0 , 0.01 );
This call integrates the system defined by harmonic_oscillator
using the RK4 method from t=0 to 10
with a step-size dt=0.01 and the initial condition
given in x
. The result,
x(t=10) is stored in x
(in-place). Each stepper defines a do_step
method which can also be used directly. So, you write down the above example
as
const double dt = 0.01; for( double t=0.0 ; t<10.0 ; t+= dt ) stepper.do_step( harmonic_oscillator , x , t , dt );
Tip | |
---|---|
If you have a C++11 enabled compiler you can easily use lambdas to create the system function :
{ runge_kutta4< state_type > stepper; integrate_const( stepper , []( const state_type &x , state_type &dxdt , double t ) { dxdt[0] = x[1]; dxdt[1] = -x[0] - gam*x[1]; } , x , 0.0 , 10.0 , 0.01 ); }
|
To improve the numerical results and additionally minimize the computational
effort, the application of a step size control is advisable. Step size
control is realized via stepper algorithms that additionally provide an
error estimation of the applied step. odeint provides a number of such
ErrorSteppers and we will show their usage
on the example of explicit_error_rk54_ck
- a 5th order Runge-Kutta method with 4th order error estimation and coefficients
introduced by Cash and Karp.
typedef runge_kutta_cash_karp54< state_type > error_stepper_type;
Given the error stepper, one still needs an instance that checks the error
and adjusts the step size accordingly. In odeint, this is done by ControlledSteppers. For the runge_kutta_cash_karp54
stepper a controlled_runge_kutta
stepper exists which can be used via
typedef controlled_runge_kutta< error_stepper_type > controlled_stepper_type; controlled_stepper_type controlled_stepper; integrate_adaptive( controlled_stepper , harmonic_oscillator , x , 0.0 , 10.0 , 0.01 );
As above, this integrates the system defined by harmonic_oscillator
,
but now using an adaptive step size method based on the Runge-Kutta Cash-Karp
54 scheme from t=0 to 10 with
an initial step size of dt=0.01 (will be adjusted)
and the initial condition given in x. The result, x(t=10),
will also be stored in x (in-place).
In the above example an error stepper is nested in a controlled stepper.
This is a nice technique; however one drawback is that one always needs
to define both steppers. One could also write the instantiation of the
controlled stepper into the call of the integrate function but a complete
knowledge of the underlying stepper types is still necessary. Another point
is, that the error tolerances for the step size control are not easily
included into the controlled stepper. Both issues can be solved by using
make_controlled
:
integrate_adaptive( make_controlled< error_stepper_type >( 1.0e-10 , 1.0e-6 ) , harmonic_oscillator , x , 0.0 , 10.0 , 0.01 );
make_controlled
can be
used with many of the steppers of odeint. The first parameter is the absolute
error tolerance eps_abs and the second is the relative
error tolerance eps_rel which is used during the integration.
The template parameter determines from which error stepper a controlled
stepper should be instantiated. An alternative syntax of make_controlled
is
integrate_adaptive( make_controlled( 1.0e-10 , 1.0e-6 , error_stepper_type() ) , harmonic_oscillator , x , 0.0 , 10.0 , 0.01 );
For the Runge-Kutta controller the error made during one step is compared
with eps_abs + eps_rel * ( ax * |x| + adxdt * dt * |dxdt| ).
If the error is smaller than this value the current step is accepted, otherwise
it is rejected and the step size is decreased. Note, that the step size
is also increased if the error gets too small compared to the rhs of the
above relation. The full instantiation of the controlled_runge_kutta
with all parameters is therefore
double abs_err = 1.0e-10 , rel_err = 1.0e-6 , a_x = 1.0 , a_dxdt = 1.0; controlled_stepper_type controlled_stepper( default_error_checker< double , range_algebra , default_operations >( abs_err , rel_err , a_x , a_dxdt ) ); integrate_adaptive( controlled_stepper , harmonic_oscillator , x , 0.0 , 10.0 , 0.01 );
When using make_controlled
the parameter ax and adxdt are
used with their standard values of 1.
In the tables below, one can find all steppers which are working with
make_controlled
and make_dense_output
which is the analog
for the dense output steppers.
Table 1.2. Generation functions make_controlled( abs_error , rel_error , stepper )
Stepper |
Result of make_controlled |
Remarks |
---|---|---|
|
|
ax=1, adxdt=1 |
|
|
ax=1, adxdt=1 |
|
|
a x=1, adxdt=1 |
|
|
- |
Table 1.3. Generation functions make_dense_output( abs_error , rel_error , stepper )
Stepper |
Result of make_dense_output |
Remarks |
---|---|---|
|
|
a x=1, adxdt=1 |
|
|
- |
When using make_controlled
or make_dense_output
one
should be aware which exact type is used and how the step size control
works.
odeint supports iterators for solving ODEs. That is, you instantiate a pair of iterators and instead of using the integrate routines with an appropriate observer you put the iterators in one of the algorithm from the C++ standard library or from Boost.Range. An example is
std::for_each( make_const_step_time_iterator_begin( stepper , harmonic_oscillator, x , 0.0 , 0.1 , 10.0 ) , make_const_step_time_iterator_end( stepper , harmonic_oscillator, x ) , []( std::pair< const state_type & , const double & > x ) { cout << x.second << " " << x.first[0] << " " << x.first[1] << "\n"; } );
The full source file for this example can be found here: harmonic_oscillator.cpp