integrate.cpp 10 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299
  1. /*
  2. [auto_generated]
  3. libs/numeric/odeint/test/integrate.cpp
  4. [begin_description]
  5. This file tests the integrate function and its variants.
  6. [end_description]
  7. Copyright 2011-2012 Mario Mulansky
  8. Copyright 2011-2012 Karsten Ahnert
  9. Distributed under the Boost Software License, Version 1.0.
  10. (See accompanying file LICENSE_1_0.txt or
  11. copy at http://www.boost.org/LICENSE_1_0.txt)
  12. */
  13. #define BOOST_TEST_MODULE odeint_integrate_functions
  14. #include <boost/type_traits.hpp>
  15. #include <vector>
  16. #include <cmath>
  17. #include <iostream>
  18. #include <boost/numeric/odeint/config.hpp>
  19. #include <boost/array.hpp>
  20. #include <boost/ref.hpp>
  21. #include <boost/iterator/counting_iterator.hpp>
  22. #include <boost/test/unit_test.hpp>
  23. #include <boost/mpl/vector.hpp>
  24. // nearly everything from odeint is used in these tests
  25. #ifndef ODEINT_INTEGRATE_ITERATOR
  26. #include <boost/numeric/odeint/integrate/integrate_const.hpp>
  27. #include <boost/numeric/odeint/integrate/integrate_adaptive.hpp>
  28. #include <boost/numeric/odeint/integrate/integrate_times.hpp>
  29. #include <boost/numeric/odeint/integrate/integrate_n_steps.hpp>
  30. #else
  31. #include <boost/numeric/odeint/iterator/integrate/integrate_const.hpp>
  32. #include <boost/numeric/odeint/iterator/integrate/integrate_adaptive.hpp>
  33. #include <boost/numeric/odeint/iterator/integrate/integrate_times.hpp>
  34. #include <boost/numeric/odeint/iterator/integrate/integrate_n_steps.hpp>
  35. #endif
  36. #include <boost/numeric/odeint/stepper/euler.hpp>
  37. #include <boost/numeric/odeint/stepper/modified_midpoint.hpp>
  38. #include <boost/numeric/odeint/stepper/runge_kutta4.hpp>
  39. #include <boost/numeric/odeint/stepper/runge_kutta_cash_karp54.hpp>
  40. #include <boost/numeric/odeint/stepper/runge_kutta_dopri5.hpp>
  41. #include <boost/numeric/odeint/stepper/runge_kutta_fehlberg78.hpp>
  42. #include <boost/numeric/odeint/stepper/controlled_runge_kutta.hpp>
  43. #include <boost/numeric/odeint/stepper/bulirsch_stoer.hpp>
  44. #include <boost/numeric/odeint/stepper/dense_output_runge_kutta.hpp>
  45. #include <boost/numeric/odeint/stepper/bulirsch_stoer_dense_out.hpp>
  46. #include <boost/numeric/odeint/stepper/adaptive_adams_bashforth_moulton.hpp>
  47. #include <boost/numeric/odeint/stepper/controlled_adams_bashforth_moulton.hpp>
  48. #include <boost/numeric/odeint/util/detail/less_with_sign.hpp>
  49. using namespace boost::unit_test;
  50. using namespace boost::numeric::odeint;
  51. namespace mpl = boost::mpl;
  52. typedef double value_type;
  53. typedef std::vector< value_type > state_type;
  54. void lorenz( const state_type &x , state_type &dxdt , const value_type t )
  55. {
  56. //const value_type sigma( 10.0 );
  57. const value_type R( 28.0 );
  58. const value_type b( value_type( 8.0 ) / value_type( 3.0 ) );
  59. // first component trivial
  60. dxdt[0] = 1.0; //sigma * ( x[1] - x[0] );
  61. dxdt[1] = R * x[0] - x[1] - x[0] * x[2];
  62. dxdt[2] = -b * x[2] + x[0] * x[1];
  63. }
  64. struct push_back_time
  65. {
  66. std::vector< double >& m_times;
  67. state_type& m_x;
  68. push_back_time( std::vector< double > &times , state_type &x )
  69. : m_times( times ) , m_x( x ) { }
  70. void operator()( const state_type &x , double t )
  71. {
  72. m_times.push_back( t );
  73. boost::numeric::odeint::copy( x , m_x );
  74. }
  75. };
  76. template< class Stepper >
  77. struct perform_integrate_const_test
  78. {
  79. void operator()( const value_type t_end , const value_type dt )
  80. {
  81. // std::cout << "Testing integrate_const with " << typeid( Stepper ).name() << std::endl;
  82. state_type x( 3 , 10.0 ) , x_end( 3 );
  83. std::vector< value_type > times;
  84. size_t steps = integrate_const( Stepper() , lorenz , x , 0.0 , t_end ,
  85. dt , push_back_time( times , x_end ) );
  86. std::cout.precision(16);
  87. std::cout << t_end << " (" << dt << "), " << steps << " , " << times.size() << " , " << 10.0+dt*steps << "=" << x_end[0] << std::endl;
  88. std::cout << static_cast<int>(floor(t_end/dt)) << " , " << t_end/dt << std::endl;
  89. BOOST_CHECK_EQUAL( static_cast<int>(times.size()) , static_cast<int>(floor(t_end/dt))+1 );
  90. for( size_t i=0 ; i<times.size() ; ++i )
  91. {
  92. //std::cout << i << " , " << times[i] << " , " << static_cast< value_type >(i)*dt << std::endl;
  93. // check if observer was called at times 0,1,2,...
  94. BOOST_CHECK_SMALL( times[i] - static_cast< value_type >(i)*dt , (i+1) * 2E-16 );
  95. }
  96. // check first, trivial, component
  97. BOOST_CHECK_SMALL( (10.0 + times[times.size()-1]) - x_end[0] , 1E-6 ); // precision of steppers: 1E-6
  98. //BOOST_CHECK_EQUAL( x[1] , x_end[1] );
  99. //BOOST_CHECK_EQUAL( x[2] , x_end[2] );
  100. }
  101. };
  102. template< class Stepper >
  103. struct perform_integrate_adaptive_test
  104. {
  105. void operator()( const value_type t_end = 10.0 , const value_type dt = 0.03 )
  106. {
  107. // std::cout << "Testing integrate_adaptive with " << typeid( Stepper ).name() << std::endl;
  108. state_type x( 3 , 10.0 ) , x_end( 3 );
  109. std::vector< value_type > times;
  110. size_t steps = integrate_adaptive( Stepper() , lorenz , x , 0.0 , t_end ,
  111. dt , push_back_time( times , x_end ) );
  112. // std::cout << t_end << " , " << steps << " , " << times.size() << " , " << dt << " , " << 10.0+t_end << "=" << x_end[0] << std::endl;
  113. BOOST_CHECK_EQUAL( times.size() , steps+1 );
  114. BOOST_CHECK_SMALL( times[0] - 0.0 , 2E-16 );
  115. BOOST_CHECK_SMALL( times[times.size()-1] - t_end , times.size() * 2E-16 );
  116. // check first, trivial, component
  117. BOOST_CHECK_SMALL( (10.0 + t_end) - x_end[0] , 1E-6 ); // precision of steppers: 1E-6
  118. // BOOST_CHECK_EQUAL( x[1] , x_end[1] );
  119. // BOOST_CHECK_EQUAL( x[2] , x_end[2] );
  120. }
  121. };
  122. template< class Stepper >
  123. struct perform_integrate_times_test
  124. {
  125. void operator()( const int n = 10 , const int dn=1 , const value_type dt = 0.03 )
  126. {
  127. // std::cout << "Testing integrate_times with " << typeid( Stepper ).name() << std::endl;
  128. state_type x( 3 ) , x_end( 3 );
  129. x[0] = x[1] = x[2] = 10.0;
  130. std::vector< double > times;
  131. std::vector< double > obs_times( abs(n) );
  132. for( int i=0 ; boost::numeric::odeint::detail::less_with_sign( static_cast<double>(i) ,
  133. static_cast<double>(obs_times.size()) ,
  134. dt ) ; i+=dn )
  135. {
  136. obs_times[i] = i;
  137. }
  138. // simple stepper
  139. integrate_times( Stepper() , lorenz , x , obs_times.begin() , obs_times.end() ,
  140. dt , push_back_time( times , x_end ) );
  141. BOOST_CHECK_EQUAL( static_cast<int>(times.size()) , abs(n) );
  142. for( size_t i=0 ; i<times.size() ; ++i )
  143. // check if observer was called at times 0,1,2,...
  144. BOOST_CHECK_EQUAL( times[i] , static_cast<double>(i) );
  145. // check first, trivial, component
  146. BOOST_CHECK_SMALL( (10.0 + 1.0*times[times.size()-1]) - x_end[0] , 1E-6 ); // precision of steppers: 1E-6
  147. // BOOST_CHECK_EQUAL( x[1] , x_end[1] );
  148. // BOOST_CHECK_EQUAL( x[2] , x_end[2] );
  149. }
  150. };
  151. template< class Stepper >
  152. struct perform_integrate_n_steps_test
  153. {
  154. void operator()( const int n = 200 , const value_type dt = 0.01 )
  155. {
  156. // std::cout << "Testing integrate_n_steps with " << typeid( Stepper ).name() << ". ";
  157. // std::cout << "dt=" << dt << std::endl;
  158. state_type x( 3 ) , x_end( 3 );
  159. x[0] = x[1] = x[2] = 10.0;
  160. std::vector< double > times;
  161. // simple stepper
  162. value_type end_time = integrate_n_steps( Stepper() , lorenz , x , 0.0 , dt , n , push_back_time( times , x_end ) );
  163. BOOST_CHECK_SMALL( end_time - n*dt , 2E-16 );
  164. BOOST_CHECK_EQUAL( static_cast<int>(times.size()) , n+1 );
  165. for( size_t i=0 ; i<times.size() ; ++i )
  166. {
  167. // check if observer was called at times 0,1,2,...
  168. BOOST_CHECK_SMALL(times[i] - static_cast< value_type >(i) * dt, 2E-16);
  169. }
  170. // check first, trivial, component
  171. BOOST_CHECK_SMALL( (10.0 + end_time) - x_end[0] , 1E-6 ); // precision of steppers: 1E-6
  172. // BOOST_CHECK_EQUAL( x[1] , x_end[1] );
  173. // BOOST_CHECK_EQUAL( x[2] , x_end[2] );
  174. }
  175. };
  176. class stepper_methods : public mpl::vector<
  177. euler< state_type > ,
  178. modified_midpoint< state_type > ,
  179. runge_kutta4< state_type > ,
  180. runge_kutta_cash_karp54< state_type > ,
  181. runge_kutta_dopri5< state_type > ,
  182. runge_kutta_fehlberg78< state_type > ,
  183. controlled_runge_kutta< runge_kutta_cash_karp54< state_type > > ,
  184. controlled_runge_kutta< runge_kutta_dopri5< state_type > > ,
  185. controlled_runge_kutta< runge_kutta_fehlberg78< state_type > > ,
  186. bulirsch_stoer< state_type > ,
  187. dense_output_runge_kutta< controlled_runge_kutta< runge_kutta_dopri5< state_type > > >,
  188. adaptive_adams_bashforth_moulton<3, state_type>,
  189. adaptive_adams_bashforth_moulton<5, state_type>,
  190. adaptive_adams_bashforth_moulton<7, state_type>,
  191. controlled_adams_bashforth_moulton<adaptive_adams_bashforth_moulton<3, state_type> >,
  192. controlled_adams_bashforth_moulton<adaptive_adams_bashforth_moulton<5, state_type> >
  193. //bulirsch_stoer_dense_out< state_type >
  194. > { };
  195. BOOST_AUTO_TEST_SUITE( integrate_test )
  196. BOOST_AUTO_TEST_CASE_TEMPLATE( integrate_const_test_case , Stepper, stepper_methods )
  197. {
  198. perform_integrate_const_test< Stepper > tester;
  199. tester( 1.005 , 0.01 );
  200. tester( 1.0 , 0.01 );
  201. tester( 1.1 , 0.01 );
  202. tester( -1.005 , -0.01 );
  203. }
  204. BOOST_AUTO_TEST_CASE_TEMPLATE( integrate_adaptive_test_case , Stepper, stepper_methods )
  205. {
  206. perform_integrate_adaptive_test< Stepper > tester;
  207. tester( 1.005 , 0.01 );
  208. tester( 1.0 , 0.01 );
  209. tester( 1.1 , 0.01 );
  210. tester( -1.005 , -0.01 );
  211. }
  212. BOOST_AUTO_TEST_CASE_TEMPLATE( integrate_times_test_case , Stepper, stepper_methods )
  213. {
  214. perform_integrate_times_test< Stepper > tester;
  215. tester();
  216. //tester( -10 , -0.01 );
  217. }
  218. BOOST_AUTO_TEST_CASE_TEMPLATE( integrate_n_steps_test_case , Stepper, stepper_methods )
  219. {
  220. if(!boost::is_same<Stepper, controlled_adams_bashforth_moulton<adaptive_adams_bashforth_moulton<3, state_type> > >::value &&
  221. !boost::is_same<Stepper, controlled_adams_bashforth_moulton<adaptive_adams_bashforth_moulton<5, state_type> > >::value)
  222. {
  223. perform_integrate_n_steps_test< Stepper > tester;
  224. tester();
  225. tester( 200 , 0.01 );
  226. tester( 200 , 0.01 );
  227. tester( 200 , 0.01 );
  228. tester( 200 , -0.01 );
  229. }
  230. }
  231. BOOST_AUTO_TEST_SUITE_END()