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