Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:42:57

0001 /*
0002  [auto_generated]
0003  boost/numeric/odeint/stepper/runge_kutta_dopri5.hpp
0004 
0005  [begin_description]
0006  Implementation of the Dormand-Prince 5(4) method. This stepper can also be used with the dense-output controlled stepper.
0007  [end_description]
0008 
0009  Copyright 2010-2013 Karsten Ahnert
0010  Copyright 2010-2013 Mario Mulansky
0011  Copyright 2012 Christoph Koke
0012 
0013  Distributed under the Boost Software License, Version 1.0.
0014  (See accompanying file LICENSE_1_0.txt or
0015  copy at http://www.boost.org/LICENSE_1_0.txt)
0016  */
0017 
0018 
0019 #ifndef BOOST_NUMERIC_ODEINT_STEPPER_RUNGE_KUTTA_DOPRI5_HPP_INCLUDED
0020 #define BOOST_NUMERIC_ODEINT_STEPPER_RUNGE_KUTTA_DOPRI5_HPP_INCLUDED
0021 
0022 
0023 #include <boost/numeric/odeint/util/bind.hpp>
0024 
0025 #include <boost/numeric/odeint/stepper/base/explicit_error_stepper_fsal_base.hpp>
0026 #include <boost/numeric/odeint/algebra/range_algebra.hpp>
0027 #include <boost/numeric/odeint/algebra/default_operations.hpp>
0028 #include <boost/numeric/odeint/algebra/algebra_dispatcher.hpp>
0029 #include <boost/numeric/odeint/algebra/operations_dispatcher.hpp>
0030 #include <boost/numeric/odeint/stepper/stepper_categories.hpp>
0031 
0032 #include <boost/numeric/odeint/util/state_wrapper.hpp>
0033 #include <boost/numeric/odeint/util/is_resizeable.hpp>
0034 #include <boost/numeric/odeint/util/resizer.hpp>
0035 #include <boost/numeric/odeint/util/same_instance.hpp>
0036 
0037 namespace boost {
0038 namespace numeric {
0039 namespace odeint {
0040 
0041 
0042 
0043 template<
0044 class State ,
0045 class Value = double ,
0046 class Deriv = State ,
0047 class Time = Value ,
0048 class Algebra = typename algebra_dispatcher< State >::algebra_type ,
0049 class Operations = typename operations_dispatcher< State >::operations_type ,
0050 class Resizer = initially_resizer
0051 >
0052 class runge_kutta_dopri5
0053 #ifndef DOXYGEN_SKIP
0054 : public explicit_error_stepper_fsal_base<
0055   runge_kutta_dopri5< State , Value , Deriv , Time , Algebra , Operations , Resizer > ,
0056   5 , 5 , 4 , State , Value , Deriv , Time , Algebra , Operations , Resizer >
0057 #else
0058 : public explicit_error_stepper_fsal_base
0059 #endif
0060 {
0061 
0062 public :
0063 
0064     #ifndef DOXYGEN_SKIP
0065     typedef explicit_error_stepper_fsal_base<
0066     runge_kutta_dopri5< State , Value , Deriv , Time , Algebra , Operations , Resizer > ,
0067     5 , 5 , 4 , State , Value , Deriv , Time , Algebra , Operations , Resizer > stepper_base_type;
0068     #else
0069     typedef explicit_error_stepper_fsal_base< runge_kutta_dopri5< ... > , ... > stepper_base_type;
0070     #endif
0071     
0072     typedef typename stepper_base_type::state_type state_type;
0073     typedef typename stepper_base_type::value_type value_type;
0074     typedef typename stepper_base_type::deriv_type deriv_type;
0075     typedef typename stepper_base_type::time_type time_type;
0076     typedef typename stepper_base_type::algebra_type algebra_type;
0077     typedef typename stepper_base_type::operations_type operations_type;
0078     typedef typename stepper_base_type::resizer_type resizer_type;
0079 
0080     #ifndef DOXYGEN_SKIP
0081     typedef typename stepper_base_type::stepper_type stepper_type;
0082     typedef typename stepper_base_type::wrapped_state_type wrapped_state_type;
0083     typedef typename stepper_base_type::wrapped_deriv_type wrapped_deriv_type;
0084     #endif // DOXYGEN_SKIP
0085 
0086 
0087     runge_kutta_dopri5( const algebra_type &algebra = algebra_type() ) : stepper_base_type( algebra )
0088     { }
0089 
0090 
0091     template< class System , class StateIn , class DerivIn , class StateOut , class DerivOut >
0092     void do_step_impl( System system , const StateIn &in , const DerivIn &dxdt_in , time_type t ,
0093             StateOut &out , DerivOut &dxdt_out , time_type dt )
0094     {
0095         const value_type a2 = static_cast<value_type> ( 1 ) / static_cast<value_type>( 5 );
0096         const value_type a3 = static_cast<value_type> ( 3 ) / static_cast<value_type> ( 10 );
0097         const value_type a4 = static_cast<value_type> ( 4 ) / static_cast<value_type> ( 5 );
0098         const value_type a5 = static_cast<value_type> ( 8 )/static_cast<value_type> ( 9 );
0099 
0100         const value_type b21 = static_cast<value_type> ( 1 ) / static_cast<value_type> ( 5 );
0101 
0102         const value_type b31 = static_cast<value_type> ( 3 ) / static_cast<value_type>( 40 );
0103         const value_type b32 = static_cast<value_type> ( 9 ) / static_cast<value_type>( 40 );
0104 
0105         const value_type b41 = static_cast<value_type> ( 44 ) / static_cast<value_type> ( 45 );
0106         const value_type b42 = static_cast<value_type> ( -56 ) / static_cast<value_type> ( 15 );
0107         const value_type b43 = static_cast<value_type> ( 32 ) / static_cast<value_type> ( 9 );
0108 
0109         const value_type b51 = static_cast<value_type> ( 19372 ) / static_cast<value_type>( 6561 );
0110         const value_type b52 = static_cast<value_type> ( -25360 ) / static_cast<value_type> ( 2187 );
0111         const value_type b53 = static_cast<value_type> ( 64448 ) / static_cast<value_type>( 6561 );
0112         const value_type b54 = static_cast<value_type> ( -212 ) / static_cast<value_type>( 729 );
0113 
0114         const value_type b61 = static_cast<value_type> ( 9017 ) / static_cast<value_type>( 3168 );
0115         const value_type b62 = static_cast<value_type> ( -355 ) / static_cast<value_type>( 33 );
0116         const value_type b63 = static_cast<value_type> ( 46732 ) / static_cast<value_type>( 5247 );
0117         const value_type b64 = static_cast<value_type> ( 49 ) / static_cast<value_type>( 176 );
0118         const value_type b65 = static_cast<value_type> ( -5103 ) / static_cast<value_type>( 18656 );
0119 
0120         const value_type c1 = static_cast<value_type> ( 35 ) / static_cast<value_type>( 384 );
0121         const value_type c3 = static_cast<value_type> ( 500 ) / static_cast<value_type>( 1113 );
0122         const value_type c4 = static_cast<value_type> ( 125 ) / static_cast<value_type>( 192 );
0123         const value_type c5 = static_cast<value_type> ( -2187 ) / static_cast<value_type>( 6784 );
0124         const value_type c6 = static_cast<value_type> ( 11 ) / static_cast<value_type>( 84 );
0125 
0126         typename odeint::unwrap_reference< System >::type &sys = system;
0127 
0128         m_k_x_tmp_resizer.adjust_size( in , detail::bind( &stepper_type::template resize_k_x_tmp_impl<StateIn> , detail::ref( *this ) , detail::_1 ) );
0129 
0130         //m_x_tmp = x + dt*b21*dxdt
0131         stepper_base_type::m_algebra.for_each3( m_x_tmp.m_v , in , dxdt_in ,
0132                 typename operations_type::template scale_sum2< value_type , time_type >( 1.0 , dt*b21 ) );
0133 
0134         sys( m_x_tmp.m_v , m_k2.m_v , t + dt*a2 );
0135         // m_x_tmp = x + dt*b31*dxdt + dt*b32*m_k2
0136         stepper_base_type::m_algebra.for_each4( m_x_tmp.m_v , in , dxdt_in , m_k2.m_v ,
0137                 typename operations_type::template scale_sum3< value_type , time_type , time_type >( 1.0 , dt*b31 , dt*b32 ));
0138 
0139         sys( m_x_tmp.m_v , m_k3.m_v , t + dt*a3 );
0140         // m_x_tmp = x + dt * (b41*dxdt + b42*m_k2 + b43*m_k3)
0141         stepper_base_type::m_algebra.for_each5( m_x_tmp.m_v , in , dxdt_in , m_k2.m_v , m_k3.m_v ,
0142                 typename operations_type::template scale_sum4< value_type , time_type , time_type , time_type >( 1.0 , dt*b41 , dt*b42 , dt*b43 ));
0143 
0144         sys( m_x_tmp.m_v, m_k4.m_v , t + dt*a4 );
0145         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 ,
0146                 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 ));
0147 
0148         sys( m_x_tmp.m_v , m_k5.m_v , t + dt*a5 );
0149         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 ,
0150                 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 ));
0151 
0152         sys( m_x_tmp.m_v , m_k6.m_v , t + dt );
0153         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 ,
0154                 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 ));
0155 
0156         // the new derivative
0157         sys( out , dxdt_out , t + dt );
0158     }
0159 
0160 
0161 
0162     template< class System , class StateIn , class DerivIn , class StateOut , class DerivOut , class Err >
0163     void do_step_impl( System system , const StateIn &in , const DerivIn &dxdt_in , time_type t ,
0164             StateOut &out , DerivOut &dxdt_out , time_type dt , Err &xerr )
0165     {
0166         const value_type c1 = static_cast<value_type> ( 35 ) / static_cast<value_type>( 384 );
0167         const value_type c3 = static_cast<value_type> ( 500 ) / static_cast<value_type>( 1113 );
0168         const value_type c4 = static_cast<value_type> ( 125 ) / static_cast<value_type>( 192 );
0169         const value_type c5 = static_cast<value_type> ( -2187 ) / static_cast<value_type>( 6784 );
0170         const value_type c6 = static_cast<value_type> ( 11 ) / static_cast<value_type>( 84 );
0171 
0172         const value_type dc1 = c1 - static_cast<value_type> ( 5179 ) / static_cast<value_type>( 57600 );
0173         const value_type dc3 = c3 - static_cast<value_type> ( 7571 ) / static_cast<value_type>( 16695 );
0174         const value_type dc4 = c4 - static_cast<value_type> ( 393 ) / static_cast<value_type>( 640 );
0175         const value_type dc5 = c5 - static_cast<value_type> ( -92097 ) / static_cast<value_type>( 339200 );
0176         const value_type dc6 = c6 - static_cast<value_type> ( 187 ) / static_cast<value_type>( 2100 );
0177         const value_type dc7 = static_cast<value_type>( -1 ) / static_cast<value_type> ( 40 );
0178 
0179         /* ToDo: copy only if &dxdt_in == &dxdt_out ? */
0180         if( same_instance( dxdt_in , dxdt_out ) )
0181         {
0182             m_dxdt_tmp_resizer.adjust_size( in , detail::bind( &stepper_type::template resize_dxdt_tmp_impl<StateIn> , detail::ref( *this ) , detail::_1 ) );
0183             boost::numeric::odeint::copy( dxdt_in , m_dxdt_tmp.m_v );
0184             do_step_impl( system , in , dxdt_in , t , out , dxdt_out , dt );
0185             //error estimate
0186             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 ,
0187                                                     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 ) );
0188 
0189         }
0190         else
0191         {
0192             do_step_impl( system , in , dxdt_in , t , out , dxdt_out , dt );
0193             //error estimate
0194             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 ,
0195                                                     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 ) );
0196         
0197         }
0198 
0199     }
0200 
0201 
0202     /*
0203      * Calculates Dense-Output for Dopri5
0204      *
0205      * See Hairer, Norsett, Wanner: Solving Ordinary Differential Equations, Nonstiff Problems. I, p.191/192
0206      *
0207      * y(t+theta) = y(t) + h * sum_i^7 b_i(theta) * k_i
0208      *
0209      * A = theta^2 * ( 3 - 2 theta )
0210      * B = theta^2 * ( theta - 1 )
0211      * C = theta^2 * ( theta - 1 )^2
0212      * D = theta   * ( theta - 1 )^2
0213      *
0214      * b_1( theta ) = A * b_1 - C * X1( theta ) + D
0215      * b_2( theta ) = 0
0216      * b_3( theta ) = A * b_3 + C * X3( theta )
0217      * b_4( theta ) = A * b_4 - C * X4( theta )
0218      * b_5( theta ) = A * b_5 + C * X5( theta )
0219      * b_6( theta ) = A * b_6 - C * X6( theta )
0220      * b_7( theta ) = B + C * X7( theta )
0221      *
0222      * An alternative Method is described in:
0223      *
0224      * www-m2.ma.tum.de/homepages/simeon/numerik3/kap3.ps
0225      */
0226     template< class StateOut , class StateIn1 , class DerivIn1 , class StateIn2 , class DerivIn2 >
0227     void calc_state( time_type t , StateOut &x ,
0228                      const StateIn1 &x_old , const DerivIn1 &deriv_old , time_type t_old ,
0229                      const StateIn2 & /* x_new */ , const DerivIn2 &deriv_new , time_type t_new ) const
0230     {
0231         const value_type b1 = static_cast<value_type> ( 35 ) / static_cast<value_type>( 384 );
0232         const value_type b3 = static_cast<value_type> ( 500 ) / static_cast<value_type>( 1113 );
0233         const value_type b4 = static_cast<value_type> ( 125 ) / static_cast<value_type>( 192 );
0234         const value_type b5 = static_cast<value_type> ( -2187 ) / static_cast<value_type>( 6784 );
0235         const value_type b6 = static_cast<value_type> ( 11 ) / static_cast<value_type>( 84 );
0236 
0237         const time_type dt = ( t_new - t_old );
0238         const value_type theta = ( t - t_old ) / dt;
0239         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 );
0240         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 );
0241         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 ) ;
0242         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 );
0243         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 );
0244         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 );
0245 
0246         const value_type theta_m_1 = theta - static_cast< value_type >( 1 );
0247         const value_type theta_sq = theta * theta;
0248         const value_type A = theta_sq * ( static_cast< value_type >( 3 ) - static_cast< value_type >( 2 ) * theta );
0249         const value_type B = theta_sq * theta_m_1;
0250         const value_type C = theta_sq * theta_m_1 * theta_m_1;
0251         const value_type D = theta * theta_m_1 * theta_m_1;
0252 
0253         const value_type b1_theta = A * b1 - C * X1 + D;
0254         const value_type b3_theta = A * b3 + C * X3;
0255         const value_type b4_theta = A * b4 - C * X4;
0256         const value_type b5_theta = A * b5 + C * X5;
0257         const value_type b6_theta = A * b6 - C * X6;
0258         const value_type b7_theta = B + C * X7;
0259 
0260         // const state_type &k1 = *m_old_deriv;
0261         // const state_type &k3 = dopri5().m_k3;
0262         // const state_type &k4 = dopri5().m_k4;
0263         // const state_type &k5 = dopri5().m_k5;
0264         // const state_type &k6 = dopri5().m_k6;
0265         // const state_type &k7 = *m_current_deriv;
0266 
0267         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 ,
0268                 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 ) );
0269     }
0270 
0271 
0272     template< class StateIn >
0273     void adjust_size( const StateIn &x )
0274     {
0275         resize_k_x_tmp_impl( x );
0276         resize_dxdt_tmp_impl( x );
0277         stepper_base_type::adjust_size( x );
0278     }
0279     
0280 
0281 private:
0282 
0283     template< class StateIn >
0284     bool resize_k_x_tmp_impl( const StateIn &x )
0285     {
0286         bool resized = false;
0287         resized |= adjust_size_by_resizeability( m_x_tmp , x , typename is_resizeable<state_type>::type() );
0288         resized |= adjust_size_by_resizeability( m_k2 , x , typename is_resizeable<deriv_type>::type() );
0289         resized |= adjust_size_by_resizeability( m_k3 , x , typename is_resizeable<deriv_type>::type() );
0290         resized |= adjust_size_by_resizeability( m_k4 , x , typename is_resizeable<deriv_type>::type() );
0291         resized |= adjust_size_by_resizeability( m_k5 , x , typename is_resizeable<deriv_type>::type() );
0292         resized |= adjust_size_by_resizeability( m_k6 , x , typename is_resizeable<deriv_type>::type() );
0293         return resized;
0294     }
0295 
0296     template< class StateIn >
0297     bool resize_dxdt_tmp_impl( const StateIn &x )
0298     {
0299         return adjust_size_by_resizeability( m_dxdt_tmp , x , typename is_resizeable<deriv_type>::type() );
0300     }
0301         
0302 
0303 
0304     wrapped_state_type m_x_tmp;
0305     wrapped_deriv_type m_k2 , m_k3 , m_k4 , m_k5 , m_k6 ;
0306     wrapped_deriv_type m_dxdt_tmp;
0307     resizer_type m_k_x_tmp_resizer;
0308     resizer_type m_dxdt_tmp_resizer;
0309 };
0310 
0311 
0312 
0313 /************* DOXYGEN ************/
0314 /**
0315  * \class runge_kutta_dopri5
0316  * \brief The Runge-Kutta Dormand-Prince 5 method.
0317  *
0318  * The Runge-Kutta Dormand-Prince 5 method is a very popular method for solving ODEs, see
0319  * <a href=""></a>.
0320  * The method is explicit and fulfills the Error Stepper concept. Step size control
0321  * is provided but continuous output is available which make this method favourable for many applications. 
0322  * 
0323  * This class derives from explicit_error_stepper_fsal_base and inherits its interface via CRTP (current recurring
0324  * template pattern). The method possesses the FSAL (first-same-as-last) property. See
0325  * explicit_error_stepper_fsal_base for more details.
0326  *
0327  * \tparam State The state type.
0328  * \tparam Value The value type.
0329  * \tparam Deriv The type representing the time derivative of the state.
0330  * \tparam Time The time representing the independent variable - the time.
0331  * \tparam Algebra The algebra type.
0332  * \tparam Operations The operations type.
0333  * \tparam Resizer The resizer policy type.
0334  */
0335 
0336 
0337     /**
0338      * \fn runge_kutta_dopri5::runge_kutta_dopri5( const algebra_type &algebra )
0339      * \brief Constructs the runge_kutta_dopri5 class. This constructor can be used as a default
0340      * constructor if the algebra has a default constructor.
0341      * \param algebra A copy of algebra is made and stored inside explicit_stepper_base.
0342      */
0343 
0344     /**
0345      * \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 )
0346      * \brief This method performs one step. The derivative `dxdt_in` of `in` at the time `t` is passed to the
0347      * method. The result is updated out-of-place, hence the input is in `in` and the output in `out`. Furthermore,
0348      * the derivative is update out-of-place, hence the input is assumed to be in `dxdt_in` and the output in
0349      * `dxdt_out`. 
0350      * Access to this step functionality is provided by explicit_error_stepper_fsal_base and 
0351      * `do_step_impl` should not be called directly.
0352      *
0353      * \param system The system function to solve, hence the r.h.s. of the ODE. It must fulfill the
0354      *               Simple System concept.
0355      * \param in The state of the ODE which should be solved. in is not modified in this method
0356      * \param dxdt_in The derivative of x at t. dxdt_in is not modified by this method
0357      * \param t The value of the time, at which the step should be performed.
0358      * \param out The result of the step is written in out.
0359      * \param dxdt_out The result of the new derivative at time t+dt.
0360      * \param dt The step size.
0361      */
0362 
0363     /**
0364      * \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 )
0365      * \brief This method performs one step. The derivative `dxdt_in` of `in` at the time `t` is passed to the
0366      * method. The result is updated out-of-place, hence the input is in `in` and the output in `out`. Furthermore,
0367      * the derivative is update out-of-place, hence the input is assumed to be in `dxdt_in` and the output in
0368      * `dxdt_out`. 
0369      * Access to this step functionality is provided by explicit_error_stepper_fsal_base and 
0370      * `do_step_impl` should not be called directly.
0371      * An estimation of the error is calculated.
0372      *
0373      * \param system The system function to solve, hence the r.h.s. of the ODE. It must fulfill the
0374      *               Simple System concept.
0375      * \param in The state of the ODE which should be solved. in is not modified in this method
0376      * \param dxdt_in The derivative of x at t. dxdt_in is not modified by this method
0377      * \param t The value of the time, at which the step should be performed.
0378      * \param out The result of the step is written in out.
0379      * \param dxdt_out The result of the new derivative at time t+dt.
0380      * \param dt The step size.
0381      * \param xerr An estimation of the error.
0382      */
0383 
0384     /**
0385      * \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
0386      * \brief This method is used for continuous output and it calculates the state `x` at a time `t` from the 
0387      * knowledge of two states `old_state` and `current_state` at time points `t_old` and `t_new`. It also uses
0388      * internal variables to calculate the result. Hence this method must be called after two successful `do_step`
0389      * calls.
0390      */
0391 
0392     /**
0393      * \fn runge_kutta_dopri5::adjust_size( const StateIn &x )
0394      * \brief Adjust the size of all temporaries in the stepper manually.
0395      * \param x A state from which the size of the temporaries to be resized is deduced.
0396      */
0397 
0398 } // odeint
0399 } // numeric
0400 } // boost
0401 
0402 
0403 #endif // BOOST_NUMERIC_ODEINT_STEPPER_RUNGE_KUTTA_DOPRI5_HPP_INCLUDED