step_size_limitation.cpp 8.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268
  1. /*
  2. [auto_generated]
  3. libs/numeric/odeint/test/step_size_limitation.cpp
  4. [begin_description]
  5. Tests the step size limitation functionality
  6. [end_description]
  7. Copyright 2015 Mario Mulansky
  8. Distributed under the Boost Software License, Version 1.0.
  9. (See accompanying file LICENSE_1_0.txt or
  10. copy at http://www.boost.org/LICENSE_1_0.txt)
  11. */
  12. #define BOOST_TEST_MODULE odeint_integrate_times
  13. #include <boost/test/unit_test.hpp>
  14. #include <utility>
  15. #include <iostream>
  16. #include <vector>
  17. #include <boost/numeric/odeint.hpp>
  18. using namespace boost::unit_test;
  19. using namespace boost::numeric::odeint;
  20. typedef double value_type;
  21. typedef std::vector< value_type > state_type;
  22. /***********************************************
  23. * first part of the tests: explicit methods
  24. ***********************************************
  25. */
  26. void damped_osc( const state_type &x , state_type &dxdt , const value_type t )
  27. {
  28. const value_type gam( 0.1);
  29. dxdt[0] = x[1];
  30. dxdt[1] = -x[0] - gam*x[1];
  31. }
  32. struct push_back_time
  33. {
  34. std::vector< double >& m_times;
  35. push_back_time( std::vector< double > &times )
  36. : m_times( times ) { }
  37. template<typename State>
  38. void operator()( const State &x , double t )
  39. {
  40. m_times.push_back( t );
  41. }
  42. };
  43. BOOST_AUTO_TEST_SUITE( step_size_limitation_test )
  44. BOOST_AUTO_TEST_CASE( test_step_adjuster )
  45. {
  46. // first use adjuster without step size limitation
  47. default_step_adjuster<double, double> step_adjuster;
  48. const double dt = 0.1;
  49. double dt_new = step_adjuster.decrease_step(dt, 1.5, 2);
  50. BOOST_CHECK(dt_new < dt*2.0/3.0);
  51. dt_new = step_adjuster.increase_step(dt, 0.8, 1);
  52. // for errors > 0.5 no increase is performed
  53. BOOST_CHECK(dt_new == dt);
  54. dt_new = step_adjuster.increase_step(dt, 0.4, 1);
  55. // smaller errors should lead to step size increase
  56. std::cout << dt_new << std::endl;
  57. BOOST_CHECK(dt_new > dt);
  58. // now test with step size limitation max_dt = 0.1
  59. default_step_adjuster<double, double>
  60. limited_adjuster(dt);
  61. dt_new = limited_adjuster.decrease_step(dt, 1.5, 2);
  62. // decreasing works as before
  63. BOOST_CHECK(dt_new < dt*2.0/3.0);
  64. dt_new = limited_adjuster.decrease_step(2*dt, 1.1, 2);
  65. // decreasing a large step size should give max_dt
  66. BOOST_CHECK(dt_new == dt);
  67. dt_new = limited_adjuster.increase_step(dt, 0.8, 1);
  68. // for errors > 0.5 no increase is performed, still valid
  69. BOOST_CHECK(dt_new == dt);
  70. dt_new = limited_adjuster.increase_step(dt, 0.4, 1);
  71. // but even for smaller errors, we should at most get 0.1
  72. BOOST_CHECK_EQUAL(dt_new, dt);
  73. dt_new = limited_adjuster.increase_step(0.9*dt, 0.1, 1);
  74. std::cout << dt_new << std::endl;
  75. // check that we don't increase beyond max_dt
  76. BOOST_CHECK(dt_new == dt);
  77. }
  78. template<class Stepper>
  79. void test_explicit_stepper(Stepper stepper, const double max_dt)
  80. {
  81. state_type x( 2 );
  82. x[0] = x[1] = 10.0;
  83. const size_t steps = 100;
  84. std::vector<double> times;
  85. integrate_adaptive(stepper, damped_osc, x, 0.0, steps*max_dt, max_dt, push_back_time(times));
  86. BOOST_CHECK_EQUAL(times.size(), steps+1);
  87. // check that dt remains at exactly max_dt
  88. for( size_t i=0 ; i<times.size() ; ++i )
  89. // check if observer was called at times 0,1,2,...
  90. BOOST_CHECK_SMALL( times[i] - static_cast<double>(i)*max_dt , 1E-15);
  91. times.clear();
  92. // this should also work when we provide some bigger initial dt
  93. integrate_adaptive(stepper, damped_osc, x, 0.0, steps*max_dt, 10*max_dt, push_back_time(times));
  94. BOOST_CHECK_EQUAL(times.size(), steps+1);
  95. // check that dt remains at exactly max_dt
  96. for( size_t i=0 ; i<times.size() ; ++i )
  97. // check if observer was called at times 0,1,2,...
  98. BOOST_CHECK_SMALL( times[i] - static_cast<double>(i)*max_dt , 1E-15);
  99. times.clear();
  100. }
  101. BOOST_AUTO_TEST_CASE( test_controlled )
  102. {
  103. const double max_dt = 0.01;
  104. test_explicit_stepper(make_controlled(1E-2, 1E-2, max_dt,
  105. runge_kutta_dopri5<state_type>()),
  106. max_dt);
  107. test_explicit_stepper(make_controlled(1E-2, 1E-2, -max_dt,
  108. runge_kutta_dopri5<state_type>()),
  109. -max_dt);
  110. test_explicit_stepper(make_controlled(1E-2, 1E-2, max_dt,
  111. runge_kutta_cash_karp54<state_type>()),
  112. max_dt);
  113. test_explicit_stepper(make_controlled(1E-2, 1E-2, -max_dt,
  114. runge_kutta_cash_karp54<state_type>()),
  115. -max_dt);
  116. test_explicit_stepper(bulirsch_stoer<state_type>(1E-2, 1E-2, 1.0, 1.0, max_dt),
  117. max_dt);
  118. test_explicit_stepper(bulirsch_stoer<state_type>(1E-2, 1E-2, 1.0, 1.0, -max_dt),
  119. -max_dt);
  120. }
  121. BOOST_AUTO_TEST_CASE( test_dense_out )
  122. {
  123. const double max_dt = 0.01;
  124. test_explicit_stepper(make_dense_output(1E-2, 1E-2, max_dt,
  125. runge_kutta_dopri5<state_type>()),
  126. max_dt);
  127. test_explicit_stepper(make_dense_output(1E-2, 1E-2, -max_dt,
  128. runge_kutta_dopri5<state_type>()),
  129. -max_dt);
  130. test_explicit_stepper(bulirsch_stoer_dense_out<state_type>(1E-2, 1E-2, 1, 1, max_dt),
  131. max_dt);
  132. test_explicit_stepper(bulirsch_stoer_dense_out<state_type>(1E-2, 1E-2, 1, 1, -max_dt),
  133. -max_dt);
  134. }
  135. /***********************************************
  136. * second part of the tests: implicit Rosenbrock
  137. ***********************************************
  138. */
  139. typedef boost::numeric::ublas::vector< value_type > vector_type;
  140. typedef boost::numeric::ublas::matrix< value_type > matrix_type;
  141. // harmonic oscillator, analytic solution x[0] = sin( t )
  142. struct osc_rhs
  143. {
  144. void operator()( const vector_type &x , vector_type &dxdt , const value_type &t ) const
  145. {
  146. dxdt( 0 ) = x( 1 );
  147. dxdt( 1 ) = -x( 0 );
  148. }
  149. };
  150. struct osc_jacobi
  151. {
  152. void operator()( const vector_type &x , matrix_type &jacobi , const value_type &t , vector_type &dfdt ) const
  153. {
  154. jacobi( 0 , 0 ) = 0;
  155. jacobi( 0 , 1 ) = 1;
  156. jacobi( 1 , 0 ) = -1;
  157. jacobi( 1 , 1 ) = 0;
  158. dfdt( 0 ) = 0.0;
  159. dfdt( 1 ) = 0.0;
  160. }
  161. };
  162. template<class Stepper>
  163. void test_rosenbrock_stepper(Stepper stepper, const double max_dt)
  164. {
  165. vector_type x( 2 );
  166. x(0) = x(1) = 10.0;
  167. const size_t steps = 100;
  168. std::vector<double> times;
  169. integrate_adaptive(stepper,
  170. std::make_pair(osc_rhs(), osc_jacobi()),
  171. x, 0.0, steps*max_dt, max_dt, push_back_time(times));
  172. BOOST_CHECK_EQUAL(times.size(), steps+1);
  173. // check that dt remains at exactly max_dt
  174. for( size_t i=0 ; i<times.size() ; ++i )
  175. // check if observer was called at times 0,1,2,...
  176. BOOST_CHECK_SMALL( times[i] - static_cast<double>(i)*max_dt , 1E-15);
  177. times.clear();
  178. // this should also work when we provide some bigger initial dt
  179. integrate_adaptive(stepper,
  180. std::make_pair(osc_rhs(), osc_jacobi()),
  181. x, 0.0, steps*max_dt, 10*max_dt, push_back_time(times));
  182. BOOST_CHECK_EQUAL(times.size(), steps+1);
  183. // check that dt remains at exactly max_dt
  184. for( size_t i=0 ; i<times.size() ; ++i )
  185. // check if observer was called at times 0,1,2,...
  186. BOOST_CHECK_SMALL( times[i] - static_cast<double>(i)*max_dt , 1E-15);
  187. times.clear();
  188. }
  189. BOOST_AUTO_TEST_CASE( test_controlled_rosenbrock )
  190. {
  191. const double max_dt = 0.01;
  192. test_rosenbrock_stepper(make_controlled(1E-2, 1E-2, max_dt, rosenbrock4<value_type>()),
  193. max_dt);
  194. test_rosenbrock_stepper(make_controlled(1E-2, 1E-2, -max_dt, rosenbrock4<value_type>()),
  195. -max_dt);
  196. }
  197. BOOST_AUTO_TEST_CASE( test_dense_out_rosenbrock )
  198. {
  199. const double max_dt = 0.01;
  200. test_rosenbrock_stepper(make_dense_output(1E-2, 1E-2, max_dt, rosenbrock4<value_type>()),
  201. max_dt);
  202. test_rosenbrock_stepper(make_dense_output(1E-2, 1E-2, -max_dt, rosenbrock4<value_type>()),
  203. -max_dt);
  204. }
  205. BOOST_AUTO_TEST_SUITE_END()