headmyshoulder / odeint-v2

odeint - solving ordinary differential equations in c++ v2
http://headmyshoulder.github.com/odeint-v2/
Other
342 stars 101 forks source link

problems with multiprecision example #109

Closed nblinov closed 11 years ago

nblinov commented 11 years ago

Hi,

I'm trying to compile lorenz_mp.cpp example that makes use of boost::multiprecision libraries but get the following error:

g++-4.8 -I./

In file included from ./boost/numeric/odeint/integrate/detail/integrate_adaptive.hpp:25:0,
                 from ./boost/numeric/odeint/integrate/integrate_adaptive.hpp:25,
                 from ./boost/numeric/odeint/integrate/integrate.hpp:24,
                 from ./boost/numeric/odeint.hpp:64,
                 from test.cpp:19:
./boost/numeric/odeint/integrate/detail/integrate_const.hpp: In instantiation of ‘size_t boost::numeric::odeint::detail::integrate_const(Stepper, System, State&, Time, Time, Time, Observer, boost::numeric::odeint::stepper_tag) [with Stepper = boost::numeric::odeint::runge_kutta4<boost::array<boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u> >, 3u>, boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u> > >; System = lorenz; State = boost::array<boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u> >, 3u>; Time = boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u> >; Observer = streaming_observer; size_t = unsigned int]’:
./boost/numeric/odeint/integrate/integrate_const.hpp:60:71:   required from ‘size_t boost::numeric::odeint::integrate_const(Stepper, System, State&, Time, Time, Time, Observer) [with Stepper = boost::numeric::odeint::runge_kutta4<boost::array<boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u> >, 3u>, boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u> > >; System = lorenz; State = boost::array<boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u> >, 3u>; Time = boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u> >; Observer = streaming_observer; size_t = unsigned int]’
test.cpp:79:40:   required from here
./boost/numeric/odeint/integrate/detail/integrate_const.hpp:54:55: error: no matching function for call to ‘less_eq_with_sign(boost::multiprecision::detail::expression<boost::multiprecision::detail::add_immediates, boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u> >, boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u> >, void, void>, boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u> >&, boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u> >&)’
     while( less_eq_with_sign( time+dt , end_time , dt ) )
                                                       ^
./boost/numeric/odeint/integrate/detail/integrate_const.hpp:54:55: note: candidate is:
In file included from ./boost/numeric/odeint/stepper/bulirsch_stoer.hpp:46:0,
                 from ./boost/numeric/odeint.hpp:39,
                 from test.cpp:19:
./boost/numeric/odeint/util/detail/less_with_sign.hpp:47:6: note: template<class T> bool boost::numeric::odeint::detail::less_eq_with_sign(T, T, T)
 bool less_eq_with_sign( T t1 , T t2 , T dt )
      ^
./boost/numeric/odeint/util/detail/less_with_sign.hpp:47:6: note:   template argument deduction/substitution failed:
In file included from ./boost/numeric/odeint/integrate/detail/integrate_adaptive.hpp:25:0,
                 from ./boost/numeric/odeint/integrate/integrate_adaptive.hpp:25,
                 from ./boost/numeric/odeint/integrate/integrate.hpp:24,
                 from ./boost/numeric/odeint.hpp:64,
                 from test.cpp:19:
./boost/numeric/odeint/integrate/detail/integrate_const.hpp:54:55: note:   deduced conflicting types for parameter ‘T’ (‘boost::multiprecision::detail::expression<boost::multiprecision::detail::add_immediates, boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u> >, boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u> >, void, void>’ and ‘boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u> >’)
     while( less_eq_with_sign( time+dt , end_time , dt ) )

I'm using odeint from the website (not the boost distro). The error is more or less the same for different gcc versions. If I change the multiprecision typedef to a double or long double everything works fine. Also, boost multiprecision examples (without odeint) compile and run fine. Any ideas of what might be going on? Thanks a lot for your help!

mariomulansky commented 11 years ago

Uhm looks like a bug i might have introduced lately, i'll have a look later.thanks for reporting.

mariomulansky commented 11 years ago

yes, my bad - and the mp examples where not built automatically so i didn't realize that earlier. I've just committed a fix, let me know if you have further problems.

nblinov commented 11 years ago

Thank you Mario! The boost multiprecision example works now (as well as the gmp-based one)!

nblinov commented 11 years ago

Hi,

in the same example, I replaced the integrate_const integrator with a plain integrate, i.e.

integrate(lorenz() , x , value_type( 0.0 ) , value_type( 10.0 ) , value_type( value_type( 1.0 ) / value_type( 10.0 ) ) ,
            streaming_observer( cout ) );

Again, if the value_type = double, there are no problems in compilation and running but if i use boost multiprecision, i get the following two errors:

In file included from ./boost/numeric/odeint/integrate/integrate_adaptive.hpp:25,
                 from ./boost/numeric/odeint/integrate/integrate.hpp:24,
                 from ./boost/numeric/odeint.hpp:64,
                 from test.cpp:19:
./boost/numeric/odeint/integrate/detail/integrate_adaptive.hpp: In function ‘size_t boost::numeric::odeint::detail::integrate_adaptive(Stepper, System, State&, Time&, Time, Time&, Observer, boost::numeric::odeint::controlled_stepper_tag) [with Stepper = boost::numeric::odeint::controlled_runge_kutta<boost::numeric::odeint::runge_kutta_dopri5<boost::array<boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u, int, void>, (boost::multiprecision::expression_template_option)1u>, 3u>, double, boost::array<boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u, int, void>, (boost::multiprecision::expression_template_option)1u>, 3u>, double, boost::numeric::odeint::array_algebra, boost::numeric::odeint::default_operations, boost::numeric::odeint::initially_resizer>, boost::numeric::odeint::default_error_checker<double, boost::numeric::odeint::array_algebra, boost::numeric::odeint::default_operations>, boost::numeric::odeint::initially_resizer, boost::numeric::odeint::explicit_error_stepper_fsal_tag>, System = lorenz, State = boost::array<boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u, int, void>, (boost::multiprecision::expression_template_option)1u>, 3u>, Time = boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u, int, void>, (boost::multiprecision::expression_template_option)1u>, Observer = streaming_observer]’:
./boost/numeric/odeint/integrate/integrate_adaptive.hpp:45:   instantiated from ‘size_t boost::numeric::odeint::integrate_adaptive(Stepper, System, State&, Time, Time, Time, Observer) [with Stepper = boost::numeric::odeint::controlled_runge_kutta<boost::numeric::odeint::runge_kutta_dopri5<boost::array<boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u, int, void>, (boost::multiprecision::expression_template_option)1u>, 3u>, double, boost::array<boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u, int, void>, (boost::multiprecision::expression_template_option)1u>, 3u>, double, boost::numeric::odeint::array_algebra, boost::numeric::odeint::default_operations, boost::numeric::odeint::initially_resizer>, boost::numeric::odeint::default_error_checker<double, boost::numeric::odeint::array_algebra, boost::numeric::odeint::default_operations>, boost::numeric::odeint::initially_resizer, boost::numeric::odeint::explicit_error_stepper_fsal_tag>, System = lorenz, State = boost::array<boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u, int, void>, (boost::multiprecision::expression_template_option)1u>, 3u>, Time = boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u, int, void>, (boost::multiprecision::expression_template_option)1u>, Observer = streaming_observer]’
./boost/numeric/odeint/integrate/integrate.hpp:42:   instantiated from ‘size_t boost::numeric::odeint::integrate(System, State&, Time, Time, Time, Observer) [with System = lorenz, State = state_type, Time = boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u, int, void>, (boost::multiprecision::expression_template_option)1u>, Observer = streaming_observer]’
test.cpp:78:   instantiated from here
./boost/numeric/odeint/integrate/detail/integrate_adaptive.hpp:92: error: no matching function for call to ‘less_with_sign(boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u, int, void>, (boost::multiprecision::expression_template_option)1u>&, boost::multiprecision::detail::expression<boost::multiprecision::detail::add_immediates, boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u, int, void>, (boost::multiprecision::expression_template_option)1u>, boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u, int, void>, (boost::multiprecision::expression_template_option)1u>, void, void>, boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u, int, void>, (boost::multiprecision::expression_template_option)1u>&)’
./boost/numeric/odeint/integrate/integrate_adaptive.hpp:45:   instantiated from ‘size_t boost::numeric::odeint::integrate_adaptive(Stepper, System, State&, Time, Time, Time, Observer) [with Stepper = boost::numeric::odeint::controlled_runge_kutta<boost::numeric::odeint::runge_kutta_dopri5<boost::array<boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u, int, void>, (boost::multiprecision::expression_template_option)1u>, 3u>, double, boost::array<boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u, int, void>, (boost::multiprecision::expression_template_option)1u>, 3u>, double, boost::numeric::odeint::array_algebra, boost::numeric::odeint::default_operations, boost::numeric::odeint::initially_resizer>, boost::numeric::odeint::default_error_checker<double, boost::numeric::odeint::array_algebra, boost::numeric::odeint::default_operations>, boost::numeric::odeint::initially_resizer, boost::numeric::odeint::explicit_error_stepper_fsal_tag>, System = lorenz, State = boost::array<boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u, int, void>, (boost::multiprecision::expression_template_option)1u>, 3u>, Time = boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u, int, void>, (boost::multiprecision::expression_template_option)1u>, Observer = streaming_observer]’
./boost/numeric/odeint/integrate/integrate.hpp:42:   instantiated from ‘size_t boost::numeric::odeint::integrate(System, State&, Time, Time, Time, Observer) [with System = lorenz, State = state_type, Time = boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u, int, void>, (boost::multiprecision::expression_template_option)1u>, Observer = streaming_observer]’
test.cpp:78:   instantiated from here
./boost/numeric/odeint/integrate/detail/integrate_adaptive.hpp:102: error: no matching function for call to ‘boost::numeric::odeint::controlled_runge_kutta<boost::numeric::odeint::runge_kutta_dopri5<boost::array<boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u, int, void>, (boost::multiprecision::expression_template_option)1u>, 3u>, double, boost::array<boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u, int, void>, (boost::multiprecision::expression_template_option)1u>, 3u>, double, boost::numeric::odeint::array_algebra, boost::numeric::odeint::default_operations, boost::numeric::odeint::initially_resizer>, boost::numeric::odeint::default_error_checker<double, boost::numeric::odeint::array_algebra, boost::numeric::odeint::default_operations>, boost::numeric::odeint::initially_resizer, boost::numeric::odeint::explicit_error_stepper_fsal_tag>::try_step(lorenz&, boost::array<boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u, int, void>, (boost::multiprecision::expression_template_option)1u>, 3u>&, boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u, int, void>, (boost::multiprecision::expression_template_option)1u>&, boost::multiprecision::number<boost::multiprecision::backends::cpp_dec_float<50u, int, void>, (boost::multiprecision::expression_template_option)1u>&)’

I think the first error (less_with_than) can be fixed by doing a explicit typecast as in the previous fixes, but I couldn't quickly figure out what the second problem is. Thanks!

mariomulansky commented 11 years ago

yes, the default template arguments where not chosen correctly. I've fixed that and also added some unit tests for mp types. However, the unit tests dont cover the integrate routines yet, so there might still be some bugs, but have a look if your code runs now.

thanks for continuously reporting!

nblinov commented 10 years ago

thanks once again! finally I'm able to do things in a reasonable time that Mathematica or scipy couldn't :)

mariomulansky commented 10 years ago

We're glad to hear odeint is useful for you! I think you can increase performance further (2x) by using a gmp backend with boost.multiprecision. Have a look at the multiprecision docs if performance is that crucial for you.