File indexing completed on 2025-11-04 09:54:23
0001 
0002 
0003 
0004 
0005 
0006 
0007 
0008 
0009 
0010 
0011 
0012 
0013 
0014 
0015 
0016 
0017 
0018 
0019 #ifndef BOOST_NUMERIC_ODEINT_STEPPER_BULIRSCH_STOER_DENSE_OUT_HPP_INCLUDED
0020 #define BOOST_NUMERIC_ODEINT_STEPPER_BULIRSCH_STOER_DENSE_OUT_HPP_INCLUDED
0021 
0022 
0023 #include <iostream>
0024 
0025 #include <algorithm>
0026 
0027 #include <boost/config.hpp> // for min/max guidelines
0028 
0029 #include <boost/numeric/odeint/util/bind.hpp>
0030 
0031 #include <boost/math/special_functions/binomial.hpp>
0032 
0033 #include <boost/numeric/odeint/stepper/controlled_runge_kutta.hpp>
0034 #include <boost/numeric/odeint/stepper/modified_midpoint.hpp>
0035 #include <boost/numeric/odeint/stepper/controlled_step_result.hpp>
0036 #include <boost/numeric/odeint/algebra/range_algebra.hpp>
0037 #include <boost/numeric/odeint/algebra/default_operations.hpp>
0038 #include <boost/numeric/odeint/algebra/algebra_dispatcher.hpp>
0039 #include <boost/numeric/odeint/algebra/operations_dispatcher.hpp>
0040 
0041 #include <boost/numeric/odeint/util/state_wrapper.hpp>
0042 #include <boost/numeric/odeint/util/is_resizeable.hpp>
0043 #include <boost/numeric/odeint/util/resizer.hpp>
0044 #include <boost/numeric/odeint/util/unit_helper.hpp>
0045 
0046 #include <boost/numeric/odeint/integrate/max_step_checker.hpp>
0047 
0048 namespace boost {
0049 namespace numeric {
0050 namespace odeint {
0051 
0052 template<
0053     class State ,
0054     class Value = double ,
0055     class Deriv = State ,
0056     class Time = Value ,
0057     class Algebra = typename algebra_dispatcher< State >::algebra_type ,
0058     class Operations = typename operations_dispatcher< State >::operations_type ,
0059     class Resizer = initially_resizer
0060     >
0061 class bulirsch_stoer_dense_out {
0062 
0063 
0064 public:
0065 
0066     typedef State state_type;
0067     typedef Value value_type;
0068     typedef Deriv deriv_type;
0069     typedef Time time_type;
0070     typedef Algebra algebra_type;
0071     typedef Operations operations_type;
0072     typedef Resizer resizer_type;
0073     typedef dense_output_stepper_tag stepper_category;
0074 #ifndef DOXYGEN_SKIP
0075     typedef state_wrapper< state_type > wrapped_state_type;
0076     typedef state_wrapper< deriv_type > wrapped_deriv_type;
0077 
0078     typedef bulirsch_stoer_dense_out< State , Value , Deriv , Time , Algebra , Operations , Resizer > controlled_error_bs_type;
0079 
0080     typedef typename inverse_time< time_type >::type inv_time_type;
0081 
0082     typedef std::vector< value_type > value_vector;
0083     typedef std::vector< time_type > time_vector;
0084     typedef std::vector< inv_time_type > inv_time_vector;  
0085     typedef std::vector< value_vector > value_matrix;
0086     typedef std::vector< size_t > int_vector;
0087     typedef std::vector< wrapped_state_type > state_vector_type;
0088     typedef std::vector< wrapped_deriv_type > deriv_vector_type;
0089     typedef std::vector< deriv_vector_type > deriv_table_type;
0090 #endif 
0091 
0092     const static size_t m_k_max = 8;
0093 
0094 
0095 
0096     bulirsch_stoer_dense_out(
0097         value_type eps_abs = 1E-6 , value_type eps_rel = 1E-6 ,
0098         value_type factor_x = 1.0 , value_type factor_dxdt = 1.0 ,
0099         time_type max_dt = static_cast<time_type>(0) ,
0100         bool control_interpolation = false )
0101         : m_error_checker( eps_abs , eps_rel , factor_x, factor_dxdt ) ,
0102           m_max_dt(max_dt) ,
0103           m_control_interpolation( control_interpolation) ,
0104           m_last_step_rejected( false ) , m_first( true ) ,
0105           m_current_state_x1( true ) ,
0106           m_error( m_k_max ) ,
0107           m_interval_sequence( m_k_max+1 ) ,
0108           m_coeff( m_k_max+1 ) ,
0109           m_cost( m_k_max+1 ) ,
0110           m_facmin_table( m_k_max+1 ) ,
0111           m_table( m_k_max ) ,
0112           m_mp_states( m_k_max+1 ) ,
0113           m_derivs( m_k_max+1 ) ,
0114           m_diffs( 2*m_k_max+2 ) ,
0115           STEPFAC1( 0.65 ) , STEPFAC2( 0.94 ) , STEPFAC3( 0.02 ) , STEPFAC4( 4.0 ) , KFAC1( 0.8 ) , KFAC2( 0.9 )
0116     {
0117         BOOST_USING_STD_MIN();
0118         BOOST_USING_STD_MAX();
0119 
0120         for( unsigned short i = 0; i < m_k_max+1; i++ )
0121         {
0122             
0123             m_interval_sequence[i] = 2 + 4*i;  
0124             m_derivs[i].resize( m_interval_sequence[i] );
0125             if( i == 0 )
0126             {
0127                 m_cost[i] = m_interval_sequence[i];
0128             } else
0129             {
0130                 m_cost[i] = m_cost[i-1] + m_interval_sequence[i];
0131             }
0132             m_facmin_table[i] = pow BOOST_PREVENT_MACRO_SUBSTITUTION( STEPFAC3 , static_cast< value_type >(1) / static_cast< value_type >( 2*i+1 ) );
0133             m_coeff[i].resize(i);
0134             for( size_t k = 0 ; k < i ; ++k  )
0135             {
0136                 const value_type r = static_cast< value_type >( m_interval_sequence[i] ) / static_cast< value_type >( m_interval_sequence[k] );
0137                 m_coeff[i][k] = 1.0 / ( r*r - static_cast< value_type >( 1.0 ) ); 
0138             }
0139             
0140 
0141             m_current_k_opt = 4;
0142             
0143 
0144 
0145 
0146         }
0147         int num = 1;
0148         for( int i = 2*(m_k_max)+1 ; i >=0  ; i-- )
0149         {
0150             m_diffs[i].resize( num );
0151             num += (i+1)%2;
0152         }
0153     }
0154 
0155     template< class System , class StateIn , class DerivIn , class StateOut , class DerivOut >
0156     controlled_step_result try_step( System system , const StateIn &in , const DerivIn &dxdt , time_type &t , StateOut &out , DerivOut &dxdt_new , time_type &dt )
0157     {
0158         if( m_max_dt != static_cast<time_type>(0) && detail::less_with_sign(m_max_dt, dt, dt) )
0159         {
0160             
0161             
0162             dt = m_max_dt;
0163             return fail;
0164         }
0165 
0166         BOOST_USING_STD_MIN();
0167         BOOST_USING_STD_MAX();
0168         using std::pow;
0169         
0170         static const value_type val1( 1.0 );
0171 
0172         bool reject( true );
0173 
0174         time_vector h_opt( m_k_max+1 );
0175         inv_time_vector work( m_k_max+1 );
0176 
0177         m_k_final = 0;
0178         time_type new_h = dt;
0179 
0180         
0181 
0182         for( size_t k = 0 ; k <= m_current_k_opt+1 ; k++ )
0183         {
0184             m_midpoint.set_steps( m_interval_sequence[k] );
0185             if( k == 0 )
0186             {
0187                 m_midpoint.do_step( system , in , dxdt , t , out , dt , m_mp_states[k].m_v , m_derivs[k]);
0188             }
0189             else
0190             {
0191                 m_midpoint.do_step( system , in , dxdt , t , m_table[k-1].m_v , dt , m_mp_states[k].m_v , m_derivs[k] );
0192                 extrapolate( k , m_table , m_coeff , out );
0193                 
0194                 m_algebra.for_each3( m_err.m_v , out , m_table[0].m_v ,
0195                                      typename operations_type::template scale_sum2< value_type , value_type >( val1 , -val1 ) );
0196                 const value_type error = m_error_checker.error( m_algebra , in , dxdt , m_err.m_v , dt );
0197                 h_opt[k] = calc_h_opt( dt , error , k );
0198                 work[k] = static_cast<value_type>( m_cost[k] ) / h_opt[k];
0199 
0200                 m_k_final = k;
0201 
0202                 if( (k == m_current_k_opt-1) || m_first )
0203                 { 
0204                     if( error < 1.0 )
0205                     {
0206                         
0207                         reject = false;
0208                         if( (work[k] < KFAC2*work[k-1]) || (m_current_k_opt <= 2) )
0209                         {
0210                             
0211                             m_current_k_opt = min BOOST_PREVENT_MACRO_SUBSTITUTION( static_cast<int>(m_k_max)-1 , max BOOST_PREVENT_MACRO_SUBSTITUTION( 2 , static_cast<int>(k)+1 ) );
0212                             new_h = h_opt[k] * static_cast<value_type>( m_cost[k+1] ) / static_cast<value_type>( m_cost[k] );
0213                         } else {
0214                             m_current_k_opt = min BOOST_PREVENT_MACRO_SUBSTITUTION( static_cast<int>(m_k_max)-1 , max BOOST_PREVENT_MACRO_SUBSTITUTION( 2 , static_cast<int>(k) ) );
0215                             new_h = h_opt[k];
0216                         }
0217                         break;
0218                     }
0219                     else if( should_reject( error , k ) && !m_first )
0220                     {
0221                         reject = true;
0222                         new_h = h_opt[k];
0223                         break;
0224                     }
0225                 }
0226                 if( k == m_current_k_opt )
0227                 { 
0228                     if( error < 1.0 )
0229                     {
0230                         
0231                         reject = false;
0232                         if( (work[k-1] < KFAC2*work[k]) )
0233                         {
0234                             m_current_k_opt = max BOOST_PREVENT_MACRO_SUBSTITUTION( 2 , static_cast<int>(m_current_k_opt)-1 );
0235                             new_h = h_opt[m_current_k_opt];
0236                         }
0237                         else if( (work[k] < KFAC2*work[k-1]) && !m_last_step_rejected )
0238                         {
0239                             m_current_k_opt = min BOOST_PREVENT_MACRO_SUBSTITUTION( static_cast<int>(m_k_max)-1 , static_cast<int>(m_current_k_opt)+1 );
0240                             new_h = h_opt[k]*static_cast<value_type>( m_cost[m_current_k_opt] ) / static_cast<value_type>( m_cost[k] );
0241                         } else
0242                             new_h = h_opt[m_current_k_opt];
0243                         break;
0244                     }
0245                     else if( should_reject( error , k ) )
0246                     {
0247                         reject = true;
0248                         new_h = h_opt[m_current_k_opt];
0249                         break;
0250                     }
0251                 }
0252                 if( k == m_current_k_opt+1 )
0253                 { 
0254                     if( error < 1.0 )
0255                     {   
0256                         reject = false;
0257                         if( work[k-2] < KFAC2*work[k-1] )
0258                             m_current_k_opt = max BOOST_PREVENT_MACRO_SUBSTITUTION( 2 , static_cast<int>(m_current_k_opt)-1 );
0259                         if( (work[k] < KFAC2*work[m_current_k_opt]) && !m_last_step_rejected )
0260                             m_current_k_opt = min BOOST_PREVENT_MACRO_SUBSTITUTION( static_cast<int>(m_k_max)-1 , static_cast<int>(k) );
0261                         new_h = h_opt[m_current_k_opt];
0262                     } else
0263                     {
0264                         reject = true;
0265                         new_h = h_opt[m_current_k_opt];
0266                     }
0267                     break;
0268                 }
0269             }
0270         }
0271 
0272         if( !reject )
0273         {
0274 
0275             
0276             typename odeint::unwrap_reference< System >::type &sys = system;
0277             sys( out , dxdt_new , t+dt );
0278 
0279             
0280             value_type error = prepare_dense_output( m_k_final , in , dxdt , out , dxdt_new , dt );
0281 
0282             if( error > static_cast<value_type>(10) ) 
0283             {
0284                 reject = true;
0285                 new_h = dt * pow BOOST_PREVENT_MACRO_SUBSTITUTION( error , static_cast<value_type>(-1)/(2*m_k_final+2) );
0286             } else {
0287                 t += dt;
0288             }
0289         }
0290         
0291         if( !m_last_step_rejected || (new_h < dt) )
0292         {
0293             
0294             if( m_max_dt != static_cast<time_type>(0) )
0295             {
0296                 new_h = detail::min_abs(m_max_dt, new_h);
0297             }
0298             dt = new_h;
0299         }
0300 
0301         m_last_step_rejected = reject;
0302         if( reject )
0303             return fail;
0304         else
0305             return success;
0306     }
0307 
0308     template< class StateType >
0309     void initialize( const StateType &x0 , const time_type &t0 , const time_type &dt0 )
0310     {
0311         m_resizer.adjust_size(x0, [this](auto&& arg) { return this->resize_impl<StateType>(std::forward<decltype(arg)>(arg)); });
0312         boost::numeric::odeint::copy( x0 , get_current_state() );
0313         m_t = t0;
0314         m_dt = dt0;
0315         reset();
0316     }
0317 
0318 
0319     
0320 
0321 
0322     template< class System >
0323     std::pair< time_type , time_type > do_step( System system )
0324     {
0325         if( m_first )
0326         {
0327             typename odeint::unwrap_reference< System >::type &sys = system;
0328             sys( get_current_state() , get_current_deriv() , m_t );
0329         }
0330 
0331         failed_step_checker fail_checker;  
0332         controlled_step_result res = fail;
0333         m_t_last = m_t;
0334         while( res == fail )
0335         {
0336             res = try_step( system , get_current_state() , get_current_deriv() , m_t , get_old_state() , get_old_deriv() , m_dt );
0337             m_first = false;
0338             fail_checker();  
0339         }
0340         toggle_current_state();
0341         return std::make_pair( m_t_last , m_t );
0342     }
0343 
0344     
0345     template< class StateOut >
0346     void calc_state( time_type t , StateOut &x ) const
0347     {
0348         do_interpolation( t , x );
0349     }
0350 
0351     const state_type& current_state( void ) const
0352     {
0353         return get_current_state();
0354     }
0355 
0356     time_type current_time( void ) const
0357     {
0358         return m_t;
0359     }
0360 
0361     const state_type& previous_state( void ) const
0362     {
0363         return get_old_state();
0364     }
0365 
0366     time_type previous_time( void ) const
0367     {
0368         return m_t_last;
0369     }
0370 
0371     time_type current_time_step( void ) const
0372     {
0373         return m_dt;
0374     }
0375 
0376     
0377     void reset()
0378     {
0379         m_first = true;
0380         m_last_step_rejected = false;
0381     }
0382 
0383     template< class StateIn >
0384     void adjust_size( const StateIn &x )
0385     {
0386         resize_impl( x );
0387         m_midpoint.adjust_size( x );
0388     }
0389 
0390 
0391 protected:
0392 
0393     time_type m_max_dt;
0394 
0395 
0396 private:
0397 
0398     template< class StateInOut , class StateVector >
0399     void extrapolate( size_t k , StateVector &table , const value_matrix &coeff , StateInOut &xest , size_t order_start_index = 0 )
0400     
0401     {
0402         static const value_type val1( 1.0 );
0403         for( int j=k-1 ; j>0 ; --j )
0404         {
0405             m_algebra.for_each3( table[j-1].m_v , table[j].m_v , table[j-1].m_v ,
0406                                  typename operations_type::template scale_sum2< value_type , value_type >( val1 + coeff[k + order_start_index][j + order_start_index] ,
0407                                                                                                            -coeff[k + order_start_index][j + order_start_index] ) );
0408         }
0409         m_algebra.for_each3( xest , table[0].m_v , xest ,
0410                              typename operations_type::template scale_sum2< value_type , value_type >( val1 + coeff[k + order_start_index][0 + order_start_index] ,
0411                                                                                                        -coeff[k + order_start_index][0 + order_start_index]) );
0412     }
0413 
0414 
0415     template< class StateVector >
0416     void extrapolate_dense_out( size_t k , StateVector &table , const value_matrix &coeff , size_t order_start_index = 0 )
0417     
0418     {
0419         
0420         static const value_type val1( 1.0 );
0421         for( int j=k ; j>1 ; --j )
0422         {
0423             m_algebra.for_each3( table[j-1].m_v , table[j].m_v , table[j-1].m_v ,
0424                                  typename operations_type::template scale_sum2< value_type , value_type >( val1 + coeff[k + order_start_index][j + order_start_index - 1] ,
0425                                                                                                            -coeff[k + order_start_index][j + order_start_index - 1] ) );
0426         }
0427         m_algebra.for_each3( table[0].m_v , table[1].m_v , table[0].m_v ,
0428                              typename operations_type::template scale_sum2< value_type , value_type >( val1 + coeff[k + order_start_index][order_start_index] ,
0429                                                                                                        -coeff[k + order_start_index][order_start_index]) );
0430     }
0431 
0432     time_type calc_h_opt( time_type h , value_type error , size_t k ) const
0433     {
0434         BOOST_USING_STD_MIN();
0435         BOOST_USING_STD_MAX();
0436         using std::pow;
0437 
0438         value_type expo = static_cast<value_type>(1)/(m_interval_sequence[k-1]);
0439         value_type facmin = m_facmin_table[k];
0440         value_type fac;
0441         if (error == 0.0)
0442             fac = static_cast<value_type>(1)/facmin;
0443         else
0444         {
0445             fac = STEPFAC2 / pow BOOST_PREVENT_MACRO_SUBSTITUTION( error / STEPFAC1 , expo );
0446             fac = max BOOST_PREVENT_MACRO_SUBSTITUTION( static_cast<value_type>( facmin/STEPFAC4 ) , min BOOST_PREVENT_MACRO_SUBSTITUTION( static_cast<value_type>(static_cast<value_type>(1)/facmin) , fac ) );
0447         }
0448         return h*fac;
0449     }
0450 
0451     bool in_convergence_window( size_t k ) const
0452     {
0453         if( (k == m_current_k_opt-1) && !m_last_step_rejected )
0454             return true; 
0455         return ( (k == m_current_k_opt) || (k == m_current_k_opt+1) );
0456     }
0457 
0458     bool should_reject( value_type error , size_t k ) const
0459     {
0460         if( k == m_current_k_opt-1 )
0461         {
0462             const value_type d = m_interval_sequence[m_current_k_opt] * m_interval_sequence[m_current_k_opt+1] /
0463                 (m_interval_sequence[0]*m_interval_sequence[0]);
0464             
0465             return ( error > d*d );
0466         }
0467         else if( k == m_current_k_opt )
0468         {
0469             const value_type d = m_interval_sequence[m_current_k_opt+1] / m_interval_sequence[0];
0470             return ( error > d*d );
0471         } else
0472             return error > 1.0;
0473     }
0474 
0475     template< class StateIn1 , class DerivIn1 , class StateIn2 , class DerivIn2 >
0476     value_type prepare_dense_output( int k , const StateIn1 &x_start , const DerivIn1 &dxdt_start ,
0477                                      const StateIn2 &  , const DerivIn2 &  , time_type dt )  
0478     
0479     {
0480 
0481         
0482 
0483 
0484 
0485 
0486 
0487 
0488 
0489         
0490         for( int j = 0 ; j<=k ; j++ )
0491         {
0492             
0493             const value_type d = m_interval_sequence[j] / ( static_cast<value_type>(2) * dt );
0494             value_type f = 1.0; 
0495             for( int kappa = 0 ; kappa <= 2*j+1 ; ++kappa )
0496             {
0497                 calculate_finite_difference( j , kappa , f , dxdt_start );
0498                 f *= d;
0499             }
0500 
0501             if( j > 0 )
0502                 extrapolate_dense_out( j , m_mp_states , m_coeff );
0503         }
0504 
0505         time_type d = dt/2;
0506 
0507         
0508         for( int kappa = 0 ; kappa<=2*k+1 ; kappa++ )
0509         {
0510             for( int j=1 ; j<=(k-kappa/2) ; ++j )
0511                 extrapolate_dense_out( j , m_diffs[kappa] , m_coeff , kappa/2 );
0512 
0513             
0514 
0515             
0516             m_algebra.for_each1( m_diffs[kappa][0].m_v , typename operations_type::template scale< time_type >( static_cast<time_type>(d) ) );
0517 
0518             d *= dt/(2*(kappa+2));
0519         }
0520 
0521         
0522 
0523         
0524         
0525 
0526         value_type error = 0.0;
0527         if( m_control_interpolation )
0528         {
0529             boost::numeric::odeint::copy( m_diffs[2*k+1][0].m_v , m_err.m_v );
0530             error = m_error_checker.error( m_algebra , x_start , dxdt_start , m_err.m_v , dt );
0531         }
0532 
0533         return error;
0534     }
0535 
0536     template< class DerivIn >
0537     void calculate_finite_difference( size_t j , size_t kappa , value_type fac , const DerivIn &dxdt )
0538     {
0539         const int m = m_interval_sequence[j]/2-1;
0540         if( kappa == 0) 
0541         {
0542             m_algebra.for_each2( m_diffs[0][j].m_v , m_derivs[j][m].m_v ,
0543                                  typename operations_type::template scale_sum1< value_type >( fac ) );
0544         }
0545         else
0546         {
0547             
0548             const int j_diffs = j - kappa/2;
0549 
0550             m_algebra.for_each2( m_diffs[kappa][j_diffs].m_v , m_derivs[j][m+kappa].m_v ,
0551                                  typename operations_type::template scale_sum1< value_type >( fac ) );
0552             value_type sign = -1.0;
0553             int c = 1;
0554             
0555             for( int i = m+static_cast<int>(kappa)-2 ; i >= m-static_cast<int>(kappa) ; i -= 2 )
0556             {
0557                 if( i >= 0 )
0558                 {
0559                     m_algebra.for_each3( m_diffs[kappa][j_diffs].m_v , m_diffs[kappa][j_diffs].m_v , m_derivs[j][i].m_v ,
0560                                          typename operations_type::template scale_sum2< value_type , value_type >( 1.0 ,
0561                                                                                                                    sign * fac * boost::math::binomial_coefficient< value_type >( kappa , c ) ) );
0562                 }
0563                 else
0564                 {
0565                     m_algebra.for_each3( m_diffs[kappa][j_diffs].m_v , m_diffs[kappa][j_diffs].m_v , dxdt ,
0566                                          typename operations_type::template scale_sum2< value_type , value_type >( 1.0 , sign * fac ) );
0567                 }
0568                 sign *= -1;
0569                 ++c;
0570             }
0571         }
0572     }
0573 
0574     template< class StateOut >
0575     void do_interpolation( time_type t , StateOut &out ) const
0576     {
0577         
0578         
0579         const value_type theta = 2 * get_unit_value( (t - m_t_last) / (m_t - m_t_last) ) - 1;
0580         
0581         
0582 
0583         boost::numeric::odeint::copy( m_mp_states[0].m_v , out );
0584         
0585         value_type theta_pow( theta );
0586         for( size_t i=0 ; i<=2*m_k_final+1 ; ++i )
0587         {
0588             m_algebra.for_each3( out , out , m_diffs[i][0].m_v ,
0589                                  typename operations_type::template scale_sum2< value_type >( static_cast<value_type>(1) , theta_pow ) );
0590             theta_pow *= theta;
0591         }
0592     }
0593 
0594     
0595     template< class StateIn >
0596     bool resize_impl( const StateIn &x )
0597     {
0598         bool resized( false );
0599 
0600         resized |= adjust_size_by_resizeability( m_x1 , x , typename is_resizeable<state_type>::type() );
0601         resized |= adjust_size_by_resizeability( m_x2 , x , typename is_resizeable<state_type>::type() );
0602         resized |= adjust_size_by_resizeability( m_dxdt1 , x , typename is_resizeable<state_type>::type() );
0603         resized |= adjust_size_by_resizeability( m_dxdt2 , x , typename is_resizeable<state_type>::type() );
0604         resized |= adjust_size_by_resizeability( m_err , x , typename is_resizeable<state_type>::type() );
0605 
0606         for( size_t i = 0 ; i < m_k_max ; ++i )
0607             resized |= adjust_size_by_resizeability( m_table[i] , x , typename is_resizeable<state_type>::type() );
0608         for( size_t i = 0 ; i < m_k_max+1 ; ++i )
0609             resized |= adjust_size_by_resizeability( m_mp_states[i] , x , typename is_resizeable<state_type>::type() );
0610         for( size_t i = 0 ; i < m_k_max+1 ; ++i )
0611             for( size_t j = 0 ; j < m_derivs[i].size() ; ++j )
0612                 resized |= adjust_size_by_resizeability( m_derivs[i][j] , x , typename is_resizeable<deriv_type>::type() );
0613         for( size_t i = 0 ; i < 2*m_k_max+2 ; ++i )
0614             for( size_t j = 0 ; j < m_diffs[i].size() ; ++j )
0615                 resized |= adjust_size_by_resizeability( m_diffs[i][j] , x , typename is_resizeable<deriv_type>::type() );
0616 
0617         return resized;
0618     }
0619 
0620 
0621     state_type& get_current_state( void )
0622     {
0623         return m_current_state_x1 ? m_x1.m_v : m_x2.m_v ;
0624     }
0625     
0626     const state_type& get_current_state( void ) const
0627     {
0628         return m_current_state_x1 ? m_x1.m_v : m_x2.m_v ;
0629     }
0630     
0631     state_type& get_old_state( void )
0632     {
0633         return m_current_state_x1 ? m_x2.m_v : m_x1.m_v ;
0634     }
0635     
0636     const state_type& get_old_state( void ) const
0637     {
0638         return m_current_state_x1 ? m_x2.m_v : m_x1.m_v ;
0639     }
0640 
0641     deriv_type& get_current_deriv( void )
0642     {
0643         return m_current_state_x1 ? m_dxdt1.m_v : m_dxdt2.m_v ;
0644     }
0645     
0646     const deriv_type& get_current_deriv( void ) const
0647     {
0648         return m_current_state_x1 ? m_dxdt1.m_v : m_dxdt2.m_v ;
0649     }
0650     
0651     deriv_type& get_old_deriv( void )
0652     {
0653         return m_current_state_x1 ? m_dxdt2.m_v : m_dxdt1.m_v ;
0654     }
0655     
0656     const deriv_type& get_old_deriv( void ) const
0657     {
0658         return m_current_state_x1 ? m_dxdt2.m_v : m_dxdt1.m_v ;
0659     }
0660 
0661     
0662     void toggle_current_state( void )
0663     {
0664         m_current_state_x1 = ! m_current_state_x1;
0665     }
0666 
0667 
0668 
0669     default_error_checker< value_type, algebra_type , operations_type > m_error_checker;
0670     modified_midpoint_dense_out< state_type , value_type , deriv_type , time_type , algebra_type , operations_type , resizer_type > m_midpoint;
0671 
0672     bool m_control_interpolation;
0673 
0674     bool m_last_step_rejected;
0675     bool m_first;
0676 
0677     time_type m_t;
0678     time_type m_dt;
0679     time_type m_dt_last;
0680     time_type m_t_last;
0681 
0682     size_t m_current_k_opt;
0683     size_t m_k_final;
0684 
0685     algebra_type m_algebra;
0686 
0687     resizer_type m_resizer;
0688 
0689     wrapped_state_type m_x1 , m_x2;
0690     wrapped_deriv_type m_dxdt1 , m_dxdt2;
0691     wrapped_state_type m_err;
0692     bool m_current_state_x1;
0693 
0694 
0695 
0696     value_vector m_error; 
0697     int_vector m_interval_sequence; 
0698     value_matrix m_coeff;
0699     int_vector m_cost; 
0700     value_vector m_facmin_table; 
0701 
0702     state_vector_type m_table; 
0703 
0704     
0705     state_vector_type m_mp_states; 
0706     deriv_table_type m_derivs; 
0707     deriv_table_type m_diffs; 
0708 
0709     
0710 
0711     value_type STEPFAC1 , STEPFAC2 , STEPFAC3 , STEPFAC4 , KFAC1 , KFAC2;
0712 };
0713 
0714 
0715 
0716 
0717 
0718 
0719 
0720 
0721 
0722 
0723 
0724 
0725 
0726 
0727 
0728 
0729 
0730 
0731 
0732 
0733 
0734 
0735 
0736     
0737 
0738 
0739 
0740 
0741 
0742 
0743 
0744 
0745 
0746 
0747 
0748 
0749     
0750 
0751 
0752 
0753 
0754 
0755 
0756 
0757 
0758 
0759 
0760 
0761 
0762 
0763 
0764 
0765 
0766 
0767 
0768 
0769 
0770     
0771 
0772 
0773 
0774 
0775 
0776 
0777 
0778 
0779     
0780 
0781 
0782 
0783 
0784 
0785 
0786 
0787 
0788 
0789 
0790     
0791 
0792 
0793 
0794 
0795 
0796 
0797 
0798     
0799 
0800 
0801 
0802 
0803 
0804     
0805 
0806 
0807 
0808 
0809 
0810     
0811 
0812 
0813 
0814 
0815 
0816     
0817 
0818 
0819 
0820 
0821 
0822     
0823 
0824 
0825 
0826 
0827 
0828     
0829 
0830 
0831 
0832 
0833 
0834 }
0835 }
0836 }
0837 
0838 #endif