Opened 7 years ago
Last modified 5 years ago
#11535 new Bugs
standard integrate_times function causes tiny time steps
Reported by: | Owned by: | karsten | |
---|---|---|---|
Milestone: | To Be Determined | Component: | odeint |
Version: | Boost 1.58.0 | Severity: | Optimization |
Keywords: | Cc: |
Description
Hello!
While implementing a new algorithm, I discovered that the function boost::numeric::odeint::detail::integrate_times(...) produces tiny time steps due to numerical accuracy issues (the relative error dt/t < 1e-18 for type long double with 18 significant digits). After changing the function as below, the runge_kutta4 algorithm became about 15% faster for my problem at the same accuracy. But maybe I am missing the reason for the old implementation?
Here is the new code with the old stuff commented out:
/* * integrate_times for simple stepper */ template< class Stepper , class System , class State , class TimeIterator , class Time , class Observer > size_t integrate_times( Stepper stepper , System system , State &start_state , TimeIterator start_time , TimeIterator end_time , Time dt , Observer observer , stepper_tag ) { BOOST_USING_STD_MIN(); typename odeint::unwrap_reference< Observer >::type &obs = observer; size_t steps = 0; Time current_dt = dt; while( true ) { Time current_time = *start_time++; obs( start_state , current_time ); if( start_time == end_time ) break; // while( less_with_sign( current_time , *start_time , current_dt ) ) // { current_dt = min BOOST_PREVENT_MACRO_SUBSTITUTION ( dt , *start_time - current_time ); stepper.do_step( system , start_state , current_time , current_dt ); // current_time += current_dt; current_time = *start_time; //THIS IS NEW! steps++; // } } return steps; }
Change History (4)
comment:1 by , 7 years ago
comment:2 by , 7 years ago
Your new code only works for the case when dt is bigger or equal to the time differences of the given time points. The while loop is exactly there to step through the solution with dt if the time points are further apart.
However, I understand your problem and I might have an idea to fix it. Nevertheless, it seems to me that you might want to look into the other integrate functions, like integrate_n_steps. Maybe they are more suitable for your case?
comment:3 by , 7 years ago
I've filed this on Github as well: https://github.com/headmyshoulder/odeint-v2/issues/171
comment:4 by , 5 years ago
Are you certain this is really a accuracy error and not related to #13524 ? You might want to check that out ...
This is the specialization of integrate_times for solvers implementing the StepperConcept. There is no accuracy for this steppers as step size control is not performed here. The loop (which you commented out) exists to reach the given time points single steps. I don't see how it can work if it is remmoved.