integrate_times.cpp 8.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247
  1. /*
  2. [auto_generated]
  3. libs/numeric/odeint/test/integrate_times.cpp
  4. [begin_description]
  5. This file tests the integrate_times function and its variants.
  6. [end_description]
  7. Copyright 2011-2012 Karsten Ahnert
  8. Copyright 2011-2012 Mario Mulansky
  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_times
  14. #include <boost/test/unit_test.hpp>
  15. #include <utility>
  16. #include <iostream>
  17. #include <vector>
  18. #include <boost/ref.hpp>
  19. #include <boost/iterator/counting_iterator.hpp>
  20. #ifndef ODEINT_INTEGRATE_ITERATOR
  21. #include <boost/numeric/odeint/integrate/integrate_times.hpp>
  22. #include <boost/numeric/odeint/integrate/integrate_adaptive.hpp>
  23. #else
  24. #include <boost/numeric/odeint/iterator/integrate/integrate_times.hpp>
  25. #include <boost/numeric/odeint/iterator/integrate/integrate_adaptive.hpp>
  26. #endif
  27. #include <boost/numeric/odeint/stepper/runge_kutta4.hpp>
  28. #include <boost/numeric/odeint/stepper/runge_kutta_dopri5.hpp>
  29. #include <boost/numeric/odeint/stepper/runge_kutta_cash_karp54.hpp>
  30. #include <boost/numeric/odeint/stepper/controlled_runge_kutta.hpp>
  31. #include <boost/numeric/odeint/stepper/bulirsch_stoer.hpp>
  32. #include <boost/numeric/odeint/stepper/bulirsch_stoer_dense_out.hpp>
  33. #include <boost/numeric/odeint/stepper/dense_output_runge_kutta.hpp>
  34. using namespace boost::unit_test;
  35. using namespace boost::numeric::odeint;
  36. typedef double value_type;
  37. typedef std::vector< value_type > state_type;
  38. void lorenz( const state_type &x , state_type &dxdt , const value_type t )
  39. {
  40. BOOST_CHECK( t >= 0.0 );
  41. const value_type sigma( 10.0 );
  42. const value_type R( 28.0 );
  43. const value_type b( value_type( 8.0 ) / value_type( 3.0 ) );
  44. dxdt[0] = sigma * ( x[1] - x[0] );
  45. dxdt[1] = R * x[0] - x[1] - x[0] * x[2];
  46. dxdt[2] = -b * x[2] + x[0] * x[1];
  47. }
  48. struct push_back_time
  49. {
  50. std::vector< double >& m_times;
  51. push_back_time( std::vector< double > &times )
  52. : m_times( times ) { }
  53. void operator()( const state_type &x , double t )
  54. {
  55. m_times.push_back( t );
  56. }
  57. };
  58. BOOST_AUTO_TEST_SUITE( integrate_times_test )
  59. BOOST_AUTO_TEST_CASE( test_integrate_times )
  60. {
  61. state_type x( 3 );
  62. x[0] = x[1] = x[2] = 10.0;
  63. const value_type dt = 0.03;
  64. std::vector< double > times;
  65. std::cout << "test rk4 stepper" << std::endl;
  66. // simple stepper
  67. integrate_times( runge_kutta4< state_type >() , lorenz , x , boost::counting_iterator<int>(0) , boost::counting_iterator<int>(10) ,
  68. dt , push_back_time( times ) );
  69. for( int i=0 ; i<10 ; ++i )
  70. // check if observer was called at times 0,1,2,...
  71. BOOST_CHECK_EQUAL( times[i] , static_cast<double>(i) );
  72. times.clear();
  73. std::cout << "test dopri5 stepper" << std::endl;
  74. // controlled stepper
  75. integrate_times( controlled_runge_kutta< runge_kutta_dopri5< state_type > >() , lorenz , x ,
  76. boost::counting_iterator<int>(0) , boost::counting_iterator<int>(10) ,
  77. dt , push_back_time( times ) );
  78. for( int i=0 ; i<10 ; ++i )
  79. // check if observer was called at times 0,1,2,...
  80. BOOST_CHECK_EQUAL( times[i] , static_cast<double>(i) );
  81. times.clear();
  82. std::cout << "test BS stepper" << std::endl;
  83. //another controlled stepper
  84. integrate_times( bulirsch_stoer< state_type >() , lorenz , x ,
  85. boost::counting_iterator<int>(0) , boost::counting_iterator<int>(10) ,
  86. dt , push_back_time( times ) );
  87. for( int i=0 ; i<10 ; ++i )
  88. // check if observer was called at times 0,1,2,...
  89. BOOST_CHECK_EQUAL( times[i] , static_cast<double>(i) );
  90. times.clear();
  91. std::cout << "test dense_out stepper" << std::endl;
  92. // dense output stepper
  93. integrate_times( dense_output_runge_kutta< controlled_runge_kutta< runge_kutta_dopri5< state_type > > >() ,
  94. lorenz , x ,
  95. boost::counting_iterator<int>(0) , boost::counting_iterator<int>(10) ,
  96. dt , push_back_time( times ) );
  97. for( int i=0 ; i<10 ; ++i )
  98. // check if observer was called at times 0,1,2,...
  99. BOOST_CHECK_EQUAL( times[i] , static_cast<double>(i) );
  100. std::cout << "test BS_do stepper" << std::endl;
  101. integrate_times( bulirsch_stoer_dense_out< state_type >() , lorenz , x ,
  102. boost::counting_iterator<int>(0) , boost::counting_iterator<int>(10) ,
  103. dt , push_back_time( times ) );
  104. for( int i=0 ; i<10 ; ++i )
  105. // check if observer was called at times 0,1,2,...
  106. BOOST_CHECK_EQUAL( times[i] , static_cast<double>(i) );
  107. }
  108. BOOST_AUTO_TEST_CASE( test_integrate_times_ranges )
  109. {
  110. state_type x( 3 );
  111. x[0] = x[1] = x[2] = 10.0;
  112. const value_type dt = 0.03;
  113. std::vector< double > times;
  114. // simple stepper
  115. integrate_times( runge_kutta4< state_type >() , lorenz , x ,
  116. std::make_pair( boost::counting_iterator<int>(0) , boost::counting_iterator<int>(10) ) ,
  117. dt , push_back_time( times ) );
  118. for( int i=0 ; i<10 ; ++i )
  119. // check if observer was called at times 0,1,2,...
  120. BOOST_CHECK_EQUAL( times[i] , static_cast<double>(i) );
  121. times.clear();
  122. // controlled stepper
  123. integrate_times( controlled_runge_kutta< runge_kutta_dopri5< state_type > >() , lorenz , x ,
  124. std::make_pair( boost::counting_iterator<int>(0) , boost::counting_iterator<int>(10) ) ,
  125. dt , push_back_time( times ) );
  126. for( int i=0 ; i<10 ; ++i )
  127. // check if observer was called at times 0,1,2,...
  128. BOOST_CHECK_EQUAL( times[i] , static_cast<double>(i) );
  129. times.clear();
  130. //another controlled stepper
  131. integrate_times( bulirsch_stoer< state_type >() , lorenz , x ,
  132. std::make_pair( boost::counting_iterator<int>(0) , boost::counting_iterator<int>(10) ) ,
  133. dt , push_back_time( times ) );
  134. for( int i=0 ; i<10 ; ++i )
  135. // check if observer was called at times 0,1,2,...
  136. BOOST_CHECK_EQUAL( times[i] , static_cast<double>(i) );
  137. times.clear();
  138. // dense output stepper
  139. integrate_times( bulirsch_stoer_dense_out< state_type >() , lorenz , x ,
  140. std::make_pair( boost::counting_iterator<int>(0) , boost::counting_iterator<int>(10) ) ,
  141. dt , push_back_time( times ) );
  142. for( int i=0 ; i<10 ; ++i )
  143. // check if observer was called at times 0,1,2,...
  144. BOOST_CHECK_EQUAL( times[i] , static_cast<double>(i) );
  145. }
  146. BOOST_AUTO_TEST_CASE( test_integrate_times_overshoot )
  147. {
  148. state_type x( 3 );
  149. x[0] = x[1] = x[2] = 10.0;
  150. double dt = -0.1;
  151. std::vector<double> times( 10 );
  152. for( int i=0 ; i<10 ; ++i )
  153. times[i] = 1.0-i*1.0/9.0;
  154. std::cout << "test rk4 stepper" << std::endl;
  155. // simple stepper
  156. std::vector<double> obs_times;
  157. int steps = integrate_times( runge_kutta4< state_type >() , lorenz , x ,
  158. times.begin() , times.end() ,
  159. dt , push_back_time( obs_times ) );
  160. // different behavior for the iterator based integrate implementaton
  161. #ifndef ODEINT_INTEGRATE_ITERATOR
  162. BOOST_CHECK_EQUAL( steps , 18 ); // we really need 18 steps because dt and
  163. // the difference of the observation times
  164. // are so out of sync
  165. #else
  166. // iterator based implementation can only return the number of iteration steps
  167. BOOST_CHECK_EQUAL( steps , 9 );
  168. #endif
  169. for( int i=0 ; i<10 ; ++i )
  170. BOOST_CHECK_EQUAL( times[i] , obs_times[i] );
  171. std::cout << "test rk_ck stepper" << std::endl;
  172. // controlled stepper
  173. obs_times.clear();
  174. integrate_times( controlled_runge_kutta< runge_kutta_cash_karp54< state_type > >() , lorenz , x ,
  175. times.begin() , times.end() ,
  176. dt , push_back_time( obs_times ) );
  177. for( int i=0 ; i<10 ; ++i )
  178. BOOST_CHECK_EQUAL( times[i] , obs_times[i] );
  179. std::cout << "test dopri5 stepper" << std::endl;
  180. // controlled stepper
  181. obs_times.clear();
  182. integrate_times( dense_output_runge_kutta< controlled_runge_kutta< runge_kutta_dopri5< state_type > > >() , lorenz , x ,
  183. times.begin() , times.end() ,
  184. dt , push_back_time( obs_times ) );
  185. for( int i=0 ; i<10 ; ++i )
  186. BOOST_CHECK_EQUAL( times[i] , obs_times[i] );
  187. }
  188. BOOST_AUTO_TEST_SUITE_END()