runge_kutta_dopri5.hpp 21 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400
  1. /*
  2. [auto_generated]
  3. boost/numeric/odeint/stepper/runge_kutta_dopri5.hpp
  4. [begin_description]
  5. Implementation of the Dormand-Prince 5(4) method. This stepper can also be used with the dense-output controlled stepper.
  6. [end_description]
  7. Copyright 2009-2011 Karsten Ahnert
  8. Copyright 2009-2011 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. #ifndef BOOST_NUMERIC_ODEINT_STEPPER_RUNGE_KUTTA_DOPRI5_HPP_INCLUDED
  14. #define BOOST_NUMERIC_ODEINT_STEPPER_RUNGE_KUTTA_DOPRI5_HPP_INCLUDED
  15. #include <boost/numeric/odeint/util/bind.hpp>
  16. #include <boost/numeric/odeint/stepper/base/explicit_error_stepper_fsal_base.hpp>
  17. #include <boost/numeric/odeint/algebra/range_algebra.hpp>
  18. #include <boost/numeric/odeint/algebra/default_operations.hpp>
  19. #include <boost/numeric/odeint/stepper/stepper_categories.hpp>
  20. #include <boost/numeric/odeint/util/state_wrapper.hpp>
  21. #include <boost/numeric/odeint/util/is_resizeable.hpp>
  22. #include <boost/numeric/odeint/util/resizer.hpp>
  23. #include <boost/numeric/odeint/util/same_instance.hpp>
  24. namespace boost {
  25. namespace numeric {
  26. namespace odeint {
  27. template<
  28. class State ,
  29. class Value = double ,
  30. class Deriv = State ,
  31. class Time = Value ,
  32. class Algebra = range_algebra ,
  33. class Operations = default_operations ,
  34. class Resizer = initially_resizer
  35. >
  36. class runge_kutta_dopri5
  37. #ifndef DOXYGEN_SKIP
  38. : public explicit_error_stepper_fsal_base<
  39. runge_kutta_dopri5< State , Value , Deriv , Time , Algebra , Operations , Resizer > ,
  40. 5 , 5 , 4 , State , Value , Deriv , Time , Algebra , Operations , Resizer >
  41. #else
  42. : public explicit_error_stepper_fsal_base
  43. #endif
  44. {
  45. public :
  46. #ifndef DOXYGEN_SKIP
  47. typedef explicit_error_stepper_fsal_base<
  48. runge_kutta_dopri5< State , Value , Deriv , Time , Algebra , Operations , Resizer > ,
  49. 5 , 5 , 4 , State , Value , Deriv , Time , Algebra , Operations , Resizer > stepper_base_type;
  50. #else
  51. typedef explicit_error_stepper_fsal_base< runge_kutta_dopri5< ... > , ... > stepper_base_type;
  52. #endif
  53. typedef typename stepper_base_type::state_type state_type;
  54. typedef typename stepper_base_type::value_type value_type;
  55. typedef typename stepper_base_type::deriv_type deriv_type;
  56. typedef typename stepper_base_type::time_type time_type;
  57. typedef typename stepper_base_type::algebra_type algebra_type;
  58. typedef typename stepper_base_type::operations_type operations_type;
  59. typedef typename stepper_base_type::resizer_type resizer_type;
  60. #ifndef DOXYGEN_SKIP
  61. typedef typename stepper_base_type::stepper_type stepper_type;
  62. typedef typename stepper_base_type::wrapped_state_type wrapped_state_type;
  63. typedef typename stepper_base_type::wrapped_deriv_type wrapped_deriv_type;
  64. #endif // DOXYGEN_SKIP
  65. runge_kutta_dopri5( const algebra_type &algebra = algebra_type() ) : stepper_base_type( algebra )
  66. { }
  67. template< class System , class StateIn , class DerivIn , class StateOut , class DerivOut >
  68. void do_step_impl( System system , const StateIn &in , const DerivIn &dxdt_in , time_type t ,
  69. StateOut &out , DerivOut &dxdt_out , time_type dt )
  70. {
  71. const value_type a2 = static_cast<value_type> ( 1 ) / static_cast<value_type>( 5 );
  72. const value_type a3 = static_cast<value_type> ( 3 ) / static_cast<value_type> ( 10 );
  73. const value_type a4 = static_cast<value_type> ( 4 ) / static_cast<value_type> ( 5 );
  74. const value_type a5 = static_cast<value_type> ( 8 )/static_cast<value_type> ( 9 );
  75. const value_type b21 = static_cast<value_type> ( 1 ) / static_cast<value_type> ( 5 );
  76. const value_type b31 = static_cast<value_type> ( 3 ) / static_cast<value_type>( 40 );
  77. const value_type b32 = static_cast<value_type> ( 9 ) / static_cast<value_type>( 40 );
  78. const value_type b41 = static_cast<value_type> ( 44 ) / static_cast<value_type> ( 45 );
  79. const value_type b42 = static_cast<value_type> ( -56 ) / static_cast<value_type> ( 15 );
  80. const value_type b43 = static_cast<value_type> ( 32 ) / static_cast<value_type> ( 9 );
  81. const value_type b51 = static_cast<value_type> ( 19372 ) / static_cast<value_type>( 6561 );
  82. const value_type b52 = static_cast<value_type> ( -25360 ) / static_cast<value_type> ( 2187 );
  83. const value_type b53 = static_cast<value_type> ( 64448 ) / static_cast<value_type>( 6561 );
  84. const value_type b54 = static_cast<value_type> ( -212 ) / static_cast<value_type>( 729 );
  85. const value_type b61 = static_cast<value_type> ( 9017 ) / static_cast<value_type>( 3168 );
  86. const value_type b62 = static_cast<value_type> ( -355 ) / static_cast<value_type>( 33 );
  87. const value_type b63 = static_cast<value_type> ( 46732 ) / static_cast<value_type>( 5247 );
  88. const value_type b64 = static_cast<value_type> ( 49 ) / static_cast<value_type>( 176 );
  89. const value_type b65 = static_cast<value_type> ( -5103 ) / static_cast<value_type>( 18656 );
  90. const value_type c1 = static_cast<value_type> ( 35 ) / static_cast<value_type>( 384 );
  91. const value_type c3 = static_cast<value_type> ( 500 ) / static_cast<value_type>( 1113 );
  92. const value_type c4 = static_cast<value_type> ( 125 ) / static_cast<value_type>( 192 );
  93. const value_type c5 = static_cast<value_type> ( -2187 ) / static_cast<value_type>( 6784 );
  94. const value_type c6 = static_cast<value_type> ( 11 ) / static_cast<value_type>( 84 );
  95. typename odeint::unwrap_reference< System >::type &sys = system;
  96. m_k_x_tmp_resizer.adjust_size( in , detail::bind( &stepper_type::template resize_k_x_tmp_impl<StateIn> , detail::ref( *this ) , detail::_1 ) );
  97. //m_x_tmp = x + dt*b21*dxdt
  98. stepper_base_type::m_algebra.for_each3( m_x_tmp.m_v , in , dxdt_in ,
  99. typename operations_type::template scale_sum2< value_type , time_type >( 1.0 , dt*b21 ) );
  100. sys( m_x_tmp.m_v , m_k2.m_v , t + dt*a2 );
  101. // m_x_tmp = x + dt*b31*dxdt + dt*b32*m_k2
  102. stepper_base_type::m_algebra.for_each4( m_x_tmp.m_v , in , dxdt_in , m_k2.m_v ,
  103. typename operations_type::template scale_sum3< value_type , time_type , time_type >( 1.0 , dt*b31 , dt*b32 ));
  104. sys( m_x_tmp.m_v , m_k3.m_v , t + dt*a3 );
  105. // m_x_tmp = x + dt * (b41*dxdt + b42*m_k2 + b43*m_k3)
  106. stepper_base_type::m_algebra.for_each5( m_x_tmp.m_v , in , dxdt_in , m_k2.m_v , m_k3.m_v ,
  107. typename operations_type::template scale_sum4< value_type , time_type , time_type , time_type >( 1.0 , dt*b41 , dt*b42 , dt*b43 ));
  108. sys( m_x_tmp.m_v, m_k4.m_v , t + dt*a4 );
  109. stepper_base_type::m_algebra.for_each6( m_x_tmp.m_v , in , dxdt_in , m_k2.m_v , m_k3.m_v , m_k4.m_v ,
  110. typename operations_type::template scale_sum5< value_type , time_type , time_type , time_type , time_type >( 1.0 , dt*b51 , dt*b52 , dt*b53 , dt*b54 ));
  111. sys( m_x_tmp.m_v , m_k5.m_v , t + dt*a5 );
  112. stepper_base_type::m_algebra.for_each7( m_x_tmp.m_v , in , dxdt_in , m_k2.m_v , m_k3.m_v , m_k4.m_v , m_k5.m_v ,
  113. typename operations_type::template scale_sum6< value_type , time_type , time_type , time_type , time_type , time_type >( 1.0 , dt*b61 , dt*b62 , dt*b63 , dt*b64 , dt*b65 ));
  114. sys( m_x_tmp.m_v , m_k6.m_v , t + dt );
  115. stepper_base_type::m_algebra.for_each7( out , in , dxdt_in , m_k3.m_v , m_k4.m_v , m_k5.m_v , m_k6.m_v ,
  116. typename operations_type::template scale_sum6< value_type , time_type , time_type , time_type , time_type , time_type >( 1.0 , dt*c1 , dt*c3 , dt*c4 , dt*c5 , dt*c6 ));
  117. // the new derivative
  118. sys( out , dxdt_out , t + dt );
  119. }
  120. template< class System , class StateIn , class DerivIn , class StateOut , class DerivOut , class Err >
  121. void do_step_impl( System system , const StateIn &in , const DerivIn &dxdt_in , time_type t ,
  122. StateOut &out , DerivOut &dxdt_out , time_type dt , Err &xerr )
  123. {
  124. const value_type c1 = static_cast<value_type> ( 35 ) / static_cast<value_type>( 384 );
  125. const value_type c3 = static_cast<value_type> ( 500 ) / static_cast<value_type>( 1113 );
  126. const value_type c4 = static_cast<value_type> ( 125 ) / static_cast<value_type>( 192 );
  127. const value_type c5 = static_cast<value_type> ( -2187 ) / static_cast<value_type>( 6784 );
  128. const value_type c6 = static_cast<value_type> ( 11 ) / static_cast<value_type>( 84 );
  129. const value_type dc1 = c1 - static_cast<value_type> ( 5179 ) / static_cast<value_type>( 57600 );
  130. const value_type dc3 = c3 - static_cast<value_type> ( 7571 ) / static_cast<value_type>( 16695 );
  131. const value_type dc4 = c4 - static_cast<value_type> ( 393 ) / static_cast<value_type>( 640 );
  132. const value_type dc5 = c5 - static_cast<value_type> ( -92097 ) / static_cast<value_type>( 339200 );
  133. const value_type dc6 = c6 - static_cast<value_type> ( 187 ) / static_cast<value_type>( 2100 );
  134. const value_type dc7 = static_cast<value_type>( -1 ) / static_cast<value_type> ( 40 );
  135. /* ToDo: copy only if &dxdt_in == &dxdt_out ? */
  136. if( same_instance( dxdt_in , dxdt_out ) )
  137. {
  138. m_dxdt_tmp_resizer.adjust_size( in , detail::bind( &stepper_type::template resize_dxdt_tmp_impl<StateIn> , detail::ref( *this ) , detail::_1 ) );
  139. boost::numeric::odeint::copy( dxdt_in , m_dxdt_tmp.m_v );
  140. do_step_impl( system , in , dxdt_in , t , out , dxdt_out , dt );
  141. //error estimate
  142. stepper_base_type::m_algebra.for_each7( xerr , m_dxdt_tmp.m_v , m_k3.m_v , m_k4.m_v , m_k5.m_v , m_k6.m_v , dxdt_out ,
  143. typename operations_type::template scale_sum6< time_type , time_type , time_type , time_type , time_type , time_type >( dt*dc1 , dt*dc3 , dt*dc4 , dt*dc5 , dt*dc6 , dt*dc7 ) );
  144. }
  145. else
  146. {
  147. do_step_impl( system , in , dxdt_in , t , out , dxdt_out , dt );
  148. //error estimate
  149. stepper_base_type::m_algebra.for_each7( xerr , dxdt_in , m_k3.m_v , m_k4.m_v , m_k5.m_v , m_k6.m_v , dxdt_out ,
  150. typename operations_type::template scale_sum6< time_type , time_type , time_type , time_type , time_type , time_type >( dt*dc1 , dt*dc3 , dt*dc4 , dt*dc5 , dt*dc6 , dt*dc7 ) );
  151. }
  152. }
  153. /*
  154. * Calculates Dense-Output for Dopri5
  155. *
  156. * See Hairer, Norsett, Wanner: Solving Ordinary Differential Equations, Nonstiff Problems. I, p.191/192
  157. *
  158. * y(t+theta) = y(t) + h * sum_i^7 b_i(theta) * k_i
  159. *
  160. * A = theta^2 * ( 3 - 2 theta )
  161. * B = theta^2 * ( theta - 1 )
  162. * C = theta^2 * ( theta - 1 )^2
  163. * D = theta * ( theta - 1 )^2
  164. *
  165. * b_1( theta ) = A * b_1 - C * X1( theta ) + D
  166. * b_2( theta ) = 0
  167. * b_3( theta ) = A * b_3 + C * X3( theta )
  168. * b_4( theta ) = A * b_4 - C * X4( theta )
  169. * b_5( theta ) = A * b_5 + C * X5( theta )
  170. * b_6( theta ) = A * b_6 - C * X6( theta )
  171. * b_7( theta ) = B + C * X7( theta )
  172. *
  173. * An alternative Method is described in:
  174. *
  175. * www-m2.ma.tum.de/homepages/simeon/numerik3/kap3.ps
  176. */
  177. template< class StateOut , class StateIn1 , class DerivIn1 , class StateIn2 , class DerivIn2 >
  178. void calc_state( time_type t , StateOut &x ,
  179. const StateIn1 &x_old , const DerivIn1 &deriv_old , time_type t_old ,
  180. const StateIn2 & /* x_new */ , const DerivIn2 &deriv_new , time_type t_new ) const
  181. {
  182. const value_type b1 = static_cast<value_type> ( 35 ) / static_cast<value_type>( 384 );
  183. const value_type b3 = static_cast<value_type> ( 500 ) / static_cast<value_type>( 1113 );
  184. const value_type b4 = static_cast<value_type> ( 125 ) / static_cast<value_type>( 192 );
  185. const value_type b5 = static_cast<value_type> ( -2187 ) / static_cast<value_type>( 6784 );
  186. const value_type b6 = static_cast<value_type> ( 11 ) / static_cast<value_type>( 84 );
  187. const time_type dt = ( t_new - t_old );
  188. const value_type theta = ( t - t_old ) / dt;
  189. const value_type X1 = static_cast< value_type >( 5 ) * ( static_cast< value_type >( 2558722523LL ) - static_cast< value_type >( 31403016 ) * theta ) / static_cast< value_type >( 11282082432LL );
  190. const value_type X3 = static_cast< value_type >( 100 ) * ( static_cast< value_type >( 882725551 ) - static_cast< value_type >( 15701508 ) * theta ) / static_cast< value_type >( 32700410799LL );
  191. const value_type X4 = static_cast< value_type >( 25 ) * ( static_cast< value_type >( 443332067 ) - static_cast< value_type >( 31403016 ) * theta ) / static_cast< value_type >( 1880347072LL ) ;
  192. const value_type X5 = static_cast< value_type >( 32805 ) * ( static_cast< value_type >( 23143187 ) - static_cast< value_type >( 3489224 ) * theta ) / static_cast< value_type >( 199316789632LL );
  193. const value_type X6 = static_cast< value_type >( 55 ) * ( static_cast< value_type >( 29972135 ) - static_cast< value_type >( 7076736 ) * theta ) / static_cast< value_type >( 822651844 );
  194. const value_type X7 = static_cast< value_type >( 10 ) * ( static_cast< value_type >( 7414447 ) - static_cast< value_type >( 829305 ) * theta ) / static_cast< value_type >( 29380423 );
  195. const value_type theta_m_1 = theta - static_cast< value_type >( 1 );
  196. const value_type theta_sq = theta * theta;
  197. const value_type A = theta_sq * ( static_cast< value_type >( 3 ) - static_cast< value_type >( 2 ) * theta );
  198. const value_type B = theta_sq * theta_m_1;
  199. const value_type C = theta_sq * theta_m_1 * theta_m_1;
  200. const value_type D = theta * theta_m_1 * theta_m_1;
  201. const value_type b1_theta = A * b1 - C * X1 + D;
  202. const value_type b3_theta = A * b3 + C * X3;
  203. const value_type b4_theta = A * b4 - C * X4;
  204. const value_type b5_theta = A * b5 + C * X5;
  205. const value_type b6_theta = A * b6 - C * X6;
  206. const value_type b7_theta = B + C * X7;
  207. // const state_type &k1 = *m_old_deriv;
  208. // const state_type &k3 = dopri5().m_k3;
  209. // const state_type &k4 = dopri5().m_k4;
  210. // const state_type &k5 = dopri5().m_k5;
  211. // const state_type &k6 = dopri5().m_k6;
  212. // const state_type &k7 = *m_current_deriv;
  213. stepper_base_type::m_algebra.for_each8( x , x_old , deriv_old , m_k3.m_v , m_k4.m_v , m_k5.m_v , m_k6.m_v , deriv_new ,
  214. typename operations_type::template scale_sum7< value_type , time_type , time_type , time_type , time_type , time_type , time_type >( 1.0 , dt * b1_theta , dt * b3_theta , dt * b4_theta , dt * b5_theta , dt * b6_theta , dt * b7_theta ) );
  215. }
  216. template< class StateIn >
  217. void adjust_size( const StateIn &x )
  218. {
  219. resize_k_x_tmp_impl( x );
  220. resize_dxdt_tmp_impl( x );
  221. stepper_base_type::adjust_size( x );
  222. }
  223. private:
  224. template< class StateIn >
  225. bool resize_k_x_tmp_impl( const StateIn &x )
  226. {
  227. bool resized = false;
  228. resized |= adjust_size_by_resizeability( m_x_tmp , x , typename is_resizeable<state_type>::type() );
  229. resized |= adjust_size_by_resizeability( m_k2 , x , typename is_resizeable<deriv_type>::type() );
  230. resized |= adjust_size_by_resizeability( m_k3 , x , typename is_resizeable<deriv_type>::type() );
  231. resized |= adjust_size_by_resizeability( m_k4 , x , typename is_resizeable<deriv_type>::type() );
  232. resized |= adjust_size_by_resizeability( m_k5 , x , typename is_resizeable<deriv_type>::type() );
  233. resized |= adjust_size_by_resizeability( m_k6 , x , typename is_resizeable<deriv_type>::type() );
  234. return resized;
  235. }
  236. template< class StateIn >
  237. bool resize_dxdt_tmp_impl( const StateIn &x )
  238. {
  239. return adjust_size_by_resizeability( m_dxdt_tmp , x , typename is_resizeable<deriv_type>::type() );
  240. }
  241. wrapped_state_type m_x_tmp;
  242. wrapped_deriv_type m_k2 , m_k3 , m_k4 , m_k5 , m_k6 ;
  243. wrapped_deriv_type m_dxdt_tmp;
  244. resizer_type m_k_x_tmp_resizer;
  245. resizer_type m_dxdt_tmp_resizer;
  246. };
  247. /************* DOXYGEN ************/
  248. /**
  249. * \class runge_kutta_dopri5
  250. * \brief The Runge-Kutta Dormand-Prince 5 method.
  251. *
  252. * The Runge-Kutta Dormand-Prince 5 method is a very popular method for solving ODEs, see
  253. * <a href=""></a>.
  254. * The method is explicit and fulfills the Error Stepper concept. Step size control
  255. * is provided but continuous output is available which make this method favourable for many applications.
  256. *
  257. * This class derives from explicit_error_stepper_fsal_base and inherits its interface via CRTP (current recurring
  258. * template pattern). The method possesses the FSAL (first-same-as-last) property. See
  259. * explicit_error_stepper_fsal_base for more details.
  260. *
  261. * \tparam State The state type.
  262. * \tparam Value The value type.
  263. * \tparam Deriv The type representing the time derivative of the state.
  264. * \tparam Time The time representing the independent variable - the time.
  265. * \tparam Algebra The algebra type.
  266. * \tparam Operations The operations type.
  267. * \tparam Resizer The resizer policy type.
  268. */
  269. /**
  270. * \fn runge_kutta_dopri5::runge_kutta_dopri5( const algebra_type &algebra )
  271. * \brief Constructs the runge_kutta_dopri5 class. This constructor can be used as a default
  272. * constructor if the algebra has a default constructor.
  273. * \param algebra A copy of algebra is made and stored inside explicit_stepper_base.
  274. */
  275. /**
  276. * \fn runge_kutta_dopri5::do_step_impl( System system , const StateIn &in , const DerivIn &dxdt_in , time_type t , StateOut &out , DerivOut &dxdt_out , time_type dt )
  277. * \brief This method performs one step. The derivative `dxdt_in` of `in` at the time `t` is passed to the
  278. * method. The result is updated out-of-place, hence the input is in `in` and the output in `out`. Furthermore,
  279. * the derivative is update out-of-place, hence the input is assumed to be in `dxdt_in` and the output in
  280. * `dxdt_out`.
  281. * Access to this step functionality is provided by explicit_error_stepper_fsal_base and
  282. * `do_step_impl` should not be called directly.
  283. *
  284. * \param system The system function to solve, hence the r.h.s. of the ODE. It must fulfill the
  285. * Simple System concept.
  286. * \param in The state of the ODE which should be solved. in is not modified in this method
  287. * \param dxdt_in The derivative of x at t. dxdt_in is not modified by this method
  288. * \param t The value of the time, at which the step should be performed.
  289. * \param out The result of the step is written in out.
  290. * \param dxdt_out The result of the new derivative at time t+dt.
  291. * \param dt The step size.
  292. */
  293. /**
  294. * \fn runge_kutta_dopri5::do_step_impl( System system , const StateIn &in , const DerivIn &dxdt_in , time_type t , StateOut &out , DerivOut &dxdt_out , time_type dt , Err &xerr )
  295. * \brief This method performs one step. The derivative `dxdt_in` of `in` at the time `t` is passed to the
  296. * method. The result is updated out-of-place, hence the input is in `in` and the output in `out`. Furthermore,
  297. * the derivative is update out-of-place, hence the input is assumed to be in `dxdt_in` and the output in
  298. * `dxdt_out`.
  299. * Access to this step functionality is provided by explicit_error_stepper_fsal_base and
  300. * `do_step_impl` should not be called directly.
  301. * An estimation of the error is calculated.
  302. *
  303. * \param system The system function to solve, hence the r.h.s. of the ODE. It must fulfill the
  304. * Simple System concept.
  305. * \param in The state of the ODE which should be solved. in is not modified in this method
  306. * \param dxdt_in The derivative of x at t. dxdt_in is not modified by this method
  307. * \param t The value of the time, at which the step should be performed.
  308. * \param out The result of the step is written in out.
  309. * \param dxdt_out The result of the new derivative at time t+dt.
  310. * \param dt The step size.
  311. * \param xerr An estimation of the error.
  312. */
  313. /**
  314. * \fn runge_kutta_dopri5::calc_state( time_type t , StateOut &x , const StateIn1 &x_old , const DerivIn1 &deriv_old , time_type t_old , const StateIn2 & , const DerivIn2 &deriv_new , time_type t_new ) const
  315. * \brief This method is used for continuous output and it calculates the state `x` at a time `t` from the
  316. * knowledge of two states `old_state` and `current_state` at time points `t_old` and `t_new`. It also uses
  317. * internal variables to calculate the result. Hence this method must be called after two successful `do_step`
  318. * calls.
  319. */
  320. /**
  321. * \fn runge_kutta_dopri5::adjust_size( const StateIn &x )
  322. * \brief Adjust the size of all temporaries in the stepper manually.
  323. * \param x A state from which the size of the temporaries to be resized is deduced.
  324. */
  325. } // odeint
  326. } // numeric
  327. } // boost
  328. #endif // BOOST_NUMERIC_ODEINT_STEPPER_RUNGE_KUTTA_DOPRI5_HPP_INCLUDED