bulirsch_stoer_dense_out.hpp 33 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841
  1. /*
  2. [auto_generated]
  3. boost/numeric/odeint/stepper/bulirsch_stoer_dense_out.hpp
  4. [begin_description]
  5. Implementaiton of the Burlish-Stoer method with dense output
  6. [end_description]
  7. Copyright 2011-2015 Mario Mulansky
  8. Copyright 2011-2013 Karsten Ahnert
  9. Copyright 2012 Christoph Koke
  10. Distributed under the Boost Software License, Version 1.0.
  11. (See accompanying file LICENSE_1_0.txt or
  12. copy at http://www.boost.org/LICENSE_1_0.txt)
  13. */
  14. #ifndef BOOST_NUMERIC_ODEINT_STEPPER_BULIRSCH_STOER_DENSE_OUT_HPP_INCLUDED
  15. #define BOOST_NUMERIC_ODEINT_STEPPER_BULIRSCH_STOER_DENSE_OUT_HPP_INCLUDED
  16. #include <iostream>
  17. #include <algorithm>
  18. #include <boost/config.hpp> // for min/max guidelines
  19. #include <boost/numeric/odeint/util/bind.hpp>
  20. #include <boost/math/special_functions/binomial.hpp>
  21. #include <boost/numeric/odeint/stepper/controlled_runge_kutta.hpp>
  22. #include <boost/numeric/odeint/stepper/modified_midpoint.hpp>
  23. #include <boost/numeric/odeint/stepper/controlled_step_result.hpp>
  24. #include <boost/numeric/odeint/algebra/range_algebra.hpp>
  25. #include <boost/numeric/odeint/algebra/default_operations.hpp>
  26. #include <boost/numeric/odeint/algebra/algebra_dispatcher.hpp>
  27. #include <boost/numeric/odeint/algebra/operations_dispatcher.hpp>
  28. #include <boost/numeric/odeint/util/state_wrapper.hpp>
  29. #include <boost/numeric/odeint/util/is_resizeable.hpp>
  30. #include <boost/numeric/odeint/util/resizer.hpp>
  31. #include <boost/numeric/odeint/util/unit_helper.hpp>
  32. #include <boost/numeric/odeint/integrate/max_step_checker.hpp>
  33. #include <boost/type_traits.hpp>
  34. namespace boost {
  35. namespace numeric {
  36. namespace odeint {
  37. template<
  38. class State ,
  39. class Value = double ,
  40. class Deriv = State ,
  41. class Time = Value ,
  42. class Algebra = typename algebra_dispatcher< State >::algebra_type ,
  43. class Operations = typename operations_dispatcher< State >::operations_type ,
  44. class Resizer = initially_resizer
  45. >
  46. class bulirsch_stoer_dense_out {
  47. public:
  48. typedef State state_type;
  49. typedef Value value_type;
  50. typedef Deriv deriv_type;
  51. typedef Time time_type;
  52. typedef Algebra algebra_type;
  53. typedef Operations operations_type;
  54. typedef Resizer resizer_type;
  55. typedef dense_output_stepper_tag stepper_category;
  56. #ifndef DOXYGEN_SKIP
  57. typedef state_wrapper< state_type > wrapped_state_type;
  58. typedef state_wrapper< deriv_type > wrapped_deriv_type;
  59. typedef bulirsch_stoer_dense_out< State , Value , Deriv , Time , Algebra , Operations , Resizer > controlled_error_bs_type;
  60. typedef typename inverse_time< time_type >::type inv_time_type;
  61. typedef std::vector< value_type > value_vector;
  62. typedef std::vector< time_type > time_vector;
  63. typedef std::vector< inv_time_type > inv_time_vector; //should be 1/time_type for boost.units
  64. typedef std::vector< value_vector > value_matrix;
  65. typedef std::vector< size_t > int_vector;
  66. typedef std::vector< wrapped_state_type > state_vector_type;
  67. typedef std::vector< wrapped_deriv_type > deriv_vector_type;
  68. typedef std::vector< deriv_vector_type > deriv_table_type;
  69. #endif //DOXYGEN_SKIP
  70. const static size_t m_k_max = 8;
  71. bulirsch_stoer_dense_out(
  72. value_type eps_abs = 1E-6 , value_type eps_rel = 1E-6 ,
  73. value_type factor_x = 1.0 , value_type factor_dxdt = 1.0 ,
  74. time_type max_dt = static_cast<time_type>(0) ,
  75. bool control_interpolation = false )
  76. : m_error_checker( eps_abs , eps_rel , factor_x, factor_dxdt ) ,
  77. m_max_dt(max_dt) ,
  78. m_control_interpolation( control_interpolation) ,
  79. m_last_step_rejected( false ) , m_first( true ) ,
  80. m_current_state_x1( true ) ,
  81. m_error( m_k_max ) ,
  82. m_interval_sequence( m_k_max+1 ) ,
  83. m_coeff( m_k_max+1 ) ,
  84. m_cost( m_k_max+1 ) ,
  85. m_facmin_table( m_k_max+1 ) ,
  86. m_table( m_k_max ) ,
  87. m_mp_states( m_k_max+1 ) ,
  88. m_derivs( m_k_max+1 ) ,
  89. m_diffs( 2*m_k_max+2 ) ,
  90. STEPFAC1( 0.65 ) , STEPFAC2( 0.94 ) , STEPFAC3( 0.02 ) , STEPFAC4( 4.0 ) , KFAC1( 0.8 ) , KFAC2( 0.9 )
  91. {
  92. BOOST_USING_STD_MIN();
  93. BOOST_USING_STD_MAX();
  94. for( unsigned short i = 0; i < m_k_max+1; i++ )
  95. {
  96. /* only this specific sequence allows for dense output */
  97. m_interval_sequence[i] = 2 + 4*i; // 2 6 10 14 ...
  98. m_derivs[i].resize( m_interval_sequence[i] );
  99. if( i == 0 )
  100. {
  101. m_cost[i] = m_interval_sequence[i];
  102. } else
  103. {
  104. m_cost[i] = m_cost[i-1] + m_interval_sequence[i];
  105. }
  106. m_facmin_table[i] = pow BOOST_PREVENT_MACRO_SUBSTITUTION( STEPFAC3 , static_cast< value_type >(1) / static_cast< value_type >( 2*i+1 ) );
  107. m_coeff[i].resize(i);
  108. for( size_t k = 0 ; k < i ; ++k )
  109. {
  110. const value_type r = static_cast< value_type >( m_interval_sequence[i] ) / static_cast< value_type >( m_interval_sequence[k] );
  111. m_coeff[i][k] = 1.0 / ( r*r - static_cast< value_type >( 1.0 ) ); // coefficients for extrapolation
  112. }
  113. // crude estimate of optimal order
  114. m_current_k_opt = 4;
  115. /* no calculation because log10 might not exist for value_type!
  116. const value_type logfact( -log10( max BOOST_PREVENT_MACRO_SUBSTITUTION( eps_rel , static_cast< value_type >( 1.0E-12 ) ) ) * 0.6 + 0.5 );
  117. m_current_k_opt = max BOOST_PREVENT_MACRO_SUBSTITUTION( 1 , min BOOST_PREVENT_MACRO_SUBSTITUTION( static_cast<int>( m_k_max-1 ) , static_cast<int>( logfact ) ));
  118. */
  119. }
  120. int num = 1;
  121. for( int i = 2*(m_k_max)+1 ; i >=0 ; i-- )
  122. {
  123. m_diffs[i].resize( num );
  124. num += (i+1)%2;
  125. }
  126. }
  127. template< class System , class StateIn , class DerivIn , class StateOut , class DerivOut >
  128. controlled_step_result try_step( System system , const StateIn &in , const DerivIn &dxdt , time_type &t , StateOut &out , DerivOut &dxdt_new , time_type &dt )
  129. {
  130. if( m_max_dt != static_cast<time_type>(0) && detail::less_with_sign(m_max_dt, dt, dt) )
  131. {
  132. // given step size is bigger then max_dt
  133. // set limit and return fail
  134. dt = m_max_dt;
  135. return fail;
  136. }
  137. BOOST_USING_STD_MIN();
  138. BOOST_USING_STD_MAX();
  139. using std::pow;
  140. static const value_type val1( 1.0 );
  141. bool reject( true );
  142. time_vector h_opt( m_k_max+1 );
  143. inv_time_vector work( m_k_max+1 );
  144. m_k_final = 0;
  145. time_type new_h = dt;
  146. //std::cout << "t=" << t <<", dt=" << dt << ", k_opt=" << m_current_k_opt << ", first: " << m_first << std::endl;
  147. for( size_t k = 0 ; k <= m_current_k_opt+1 ; k++ )
  148. {
  149. m_midpoint.set_steps( m_interval_sequence[k] );
  150. if( k == 0 )
  151. {
  152. m_midpoint.do_step( system , in , dxdt , t , out , dt , m_mp_states[k].m_v , m_derivs[k]);
  153. }
  154. else
  155. {
  156. m_midpoint.do_step( system , in , dxdt , t , m_table[k-1].m_v , dt , m_mp_states[k].m_v , m_derivs[k] );
  157. extrapolate( k , m_table , m_coeff , out );
  158. // get error estimate
  159. m_algebra.for_each3( m_err.m_v , out , m_table[0].m_v ,
  160. typename operations_type::template scale_sum2< value_type , value_type >( val1 , -val1 ) );
  161. const value_type error = m_error_checker.error( m_algebra , in , dxdt , m_err.m_v , dt );
  162. h_opt[k] = calc_h_opt( dt , error , k );
  163. work[k] = static_cast<value_type>( m_cost[k] ) / h_opt[k];
  164. m_k_final = k;
  165. if( (k == m_current_k_opt-1) || m_first )
  166. { // convergence before k_opt ?
  167. if( error < 1.0 )
  168. {
  169. //convergence
  170. reject = false;
  171. if( (work[k] < KFAC2*work[k-1]) || (m_current_k_opt <= 2) )
  172. {
  173. // leave order as is (except we were in first round)
  174. 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 ) );
  175. new_h = h_opt[k] * static_cast<value_type>( m_cost[k+1] ) / static_cast<value_type>( m_cost[k] );
  176. } else {
  177. 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) ) );
  178. new_h = h_opt[k];
  179. }
  180. break;
  181. }
  182. else if( should_reject( error , k ) && !m_first )
  183. {
  184. reject = true;
  185. new_h = h_opt[k];
  186. break;
  187. }
  188. }
  189. if( k == m_current_k_opt )
  190. { // convergence at k_opt ?
  191. if( error < 1.0 )
  192. {
  193. //convergence
  194. reject = false;
  195. if( (work[k-1] < KFAC2*work[k]) )
  196. {
  197. m_current_k_opt = max BOOST_PREVENT_MACRO_SUBSTITUTION( 2 , static_cast<int>(m_current_k_opt)-1 );
  198. new_h = h_opt[m_current_k_opt];
  199. }
  200. else if( (work[k] < KFAC2*work[k-1]) && !m_last_step_rejected )
  201. {
  202. m_current_k_opt = min BOOST_PREVENT_MACRO_SUBSTITUTION( static_cast<int>(m_k_max)-1 , static_cast<int>(m_current_k_opt)+1 );
  203. new_h = h_opt[k]*static_cast<value_type>( m_cost[m_current_k_opt] ) / static_cast<value_type>( m_cost[k] );
  204. } else
  205. new_h = h_opt[m_current_k_opt];
  206. break;
  207. }
  208. else if( should_reject( error , k ) )
  209. {
  210. reject = true;
  211. new_h = h_opt[m_current_k_opt];
  212. break;
  213. }
  214. }
  215. if( k == m_current_k_opt+1 )
  216. { // convergence at k_opt+1 ?
  217. if( error < 1.0 )
  218. { //convergence
  219. reject = false;
  220. if( work[k-2] < KFAC2*work[k-1] )
  221. m_current_k_opt = max BOOST_PREVENT_MACRO_SUBSTITUTION( 2 , static_cast<int>(m_current_k_opt)-1 );
  222. if( (work[k] < KFAC2*work[m_current_k_opt]) && !m_last_step_rejected )
  223. m_current_k_opt = min BOOST_PREVENT_MACRO_SUBSTITUTION( static_cast<int>(m_k_max)-1 , static_cast<int>(k) );
  224. new_h = h_opt[m_current_k_opt];
  225. } else
  226. {
  227. reject = true;
  228. new_h = h_opt[m_current_k_opt];
  229. }
  230. break;
  231. }
  232. }
  233. }
  234. if( !reject )
  235. {
  236. //calculate dxdt for next step and dense output
  237. typename odeint::unwrap_reference< System >::type &sys = system;
  238. sys( out , dxdt_new , t+dt );
  239. //prepare dense output
  240. value_type error = prepare_dense_output( m_k_final , in , dxdt , out , dxdt_new , dt );
  241. if( error > static_cast<value_type>(10) ) // we are not as accurate for interpolation as for the steps
  242. {
  243. reject = true;
  244. new_h = dt * pow BOOST_PREVENT_MACRO_SUBSTITUTION( error , static_cast<value_type>(-1)/(2*m_k_final+2) );
  245. } else {
  246. t += dt;
  247. }
  248. }
  249. //set next stepsize
  250. if( !m_last_step_rejected || (new_h < dt) )
  251. {
  252. // limit step size
  253. if( m_max_dt != static_cast<time_type>(0) )
  254. {
  255. new_h = detail::min_abs(m_max_dt, new_h);
  256. }
  257. dt = new_h;
  258. }
  259. m_last_step_rejected = reject;
  260. if( reject )
  261. return fail;
  262. else
  263. return success;
  264. }
  265. template< class StateType >
  266. void initialize( const StateType &x0 , const time_type &t0 , const time_type &dt0 )
  267. {
  268. m_resizer.adjust_size( x0 , detail::bind( &controlled_error_bs_type::template resize_impl< StateType > , detail::ref( *this ) , detail::_1 ) );
  269. boost::numeric::odeint::copy( x0 , get_current_state() );
  270. m_t = t0;
  271. m_dt = dt0;
  272. reset();
  273. }
  274. /* =======================================================
  275. * the actual step method that should be called from outside (maybe make try_step private?)
  276. */
  277. template< class System >
  278. std::pair< time_type , time_type > do_step( System system )
  279. {
  280. if( m_first )
  281. {
  282. typename odeint::unwrap_reference< System >::type &sys = system;
  283. sys( get_current_state() , get_current_deriv() , m_t );
  284. }
  285. failed_step_checker fail_checker; // to throw a runtime_error if step size adjustment fails
  286. controlled_step_result res = fail;
  287. m_t_last = m_t;
  288. while( res == fail )
  289. {
  290. res = try_step( system , get_current_state() , get_current_deriv() , m_t , get_old_state() , get_old_deriv() , m_dt );
  291. m_first = false;
  292. fail_checker(); // check for overflow of failed steps
  293. }
  294. toggle_current_state();
  295. return std::make_pair( m_t_last , m_t );
  296. }
  297. /* performs the interpolation from a calculated step */
  298. template< class StateOut >
  299. void calc_state( time_type t , StateOut &x ) const
  300. {
  301. do_interpolation( t , x );
  302. }
  303. const state_type& current_state( void ) const
  304. {
  305. return get_current_state();
  306. }
  307. time_type current_time( void ) const
  308. {
  309. return m_t;
  310. }
  311. const state_type& previous_state( void ) const
  312. {
  313. return get_old_state();
  314. }
  315. time_type previous_time( void ) const
  316. {
  317. return m_t_last;
  318. }
  319. time_type current_time_step( void ) const
  320. {
  321. return m_dt;
  322. }
  323. /** \brief Resets the internal state of the stepper. */
  324. void reset()
  325. {
  326. m_first = true;
  327. m_last_step_rejected = false;
  328. }
  329. template< class StateIn >
  330. void adjust_size( const StateIn &x )
  331. {
  332. resize_impl( x );
  333. m_midpoint.adjust_size( x );
  334. }
  335. protected:
  336. time_type m_max_dt;
  337. private:
  338. template< class StateInOut , class StateVector >
  339. void extrapolate( size_t k , StateVector &table , const value_matrix &coeff , StateInOut &xest , size_t order_start_index = 0 )
  340. //polynomial extrapolation, see http://www.nr.com/webnotes/nr3web21.pdf
  341. {
  342. static const value_type val1( 1.0 );
  343. for( int j=k-1 ; j>0 ; --j )
  344. {
  345. m_algebra.for_each3( table[j-1].m_v , table[j].m_v , table[j-1].m_v ,
  346. typename operations_type::template scale_sum2< value_type , value_type >( val1 + coeff[k + order_start_index][j + order_start_index] ,
  347. -coeff[k + order_start_index][j + order_start_index] ) );
  348. }
  349. m_algebra.for_each3( xest , table[0].m_v , xest ,
  350. typename operations_type::template scale_sum2< value_type , value_type >( val1 + coeff[k + order_start_index][0 + order_start_index] ,
  351. -coeff[k + order_start_index][0 + order_start_index]) );
  352. }
  353. template< class StateVector >
  354. void extrapolate_dense_out( size_t k , StateVector &table , const value_matrix &coeff , size_t order_start_index = 0 )
  355. //polynomial extrapolation, see http://www.nr.com/webnotes/nr3web21.pdf
  356. {
  357. // result is written into table[0]
  358. static const value_type val1( 1.0 );
  359. for( int j=k ; j>1 ; --j )
  360. {
  361. m_algebra.for_each3( table[j-1].m_v , table[j].m_v , table[j-1].m_v ,
  362. typename operations_type::template scale_sum2< value_type , value_type >( val1 + coeff[k + order_start_index][j + order_start_index - 1] ,
  363. -coeff[k + order_start_index][j + order_start_index - 1] ) );
  364. }
  365. m_algebra.for_each3( table[0].m_v , table[1].m_v , table[0].m_v ,
  366. typename operations_type::template scale_sum2< value_type , value_type >( val1 + coeff[k + order_start_index][order_start_index] ,
  367. -coeff[k + order_start_index][order_start_index]) );
  368. }
  369. time_type calc_h_opt( time_type h , value_type error , size_t k ) const
  370. {
  371. BOOST_USING_STD_MIN();
  372. BOOST_USING_STD_MAX();
  373. using std::pow;
  374. value_type expo = static_cast<value_type>(1)/(m_interval_sequence[k-1]);
  375. value_type facmin = m_facmin_table[k];
  376. value_type fac;
  377. if (error == 0.0)
  378. fac = static_cast<value_type>(1)/facmin;
  379. else
  380. {
  381. fac = STEPFAC2 / pow BOOST_PREVENT_MACRO_SUBSTITUTION( error / STEPFAC1 , expo );
  382. 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 ) );
  383. }
  384. return h*fac;
  385. }
  386. bool in_convergence_window( size_t k ) const
  387. {
  388. if( (k == m_current_k_opt-1) && !m_last_step_rejected )
  389. return true; // decrease order only if last step was not rejected
  390. return ( (k == m_current_k_opt) || (k == m_current_k_opt+1) );
  391. }
  392. bool should_reject( value_type error , size_t k ) const
  393. {
  394. if( k == m_current_k_opt-1 )
  395. {
  396. const value_type d = m_interval_sequence[m_current_k_opt] * m_interval_sequence[m_current_k_opt+1] /
  397. (m_interval_sequence[0]*m_interval_sequence[0]);
  398. //step will fail, criterion 17.3.17 in NR
  399. return ( error > d*d );
  400. }
  401. else if( k == m_current_k_opt )
  402. {
  403. const value_type d = m_interval_sequence[m_current_k_opt+1] / m_interval_sequence[0];
  404. return ( error > d*d );
  405. } else
  406. return error > 1.0;
  407. }
  408. template< class StateIn1 , class DerivIn1 , class StateIn2 , class DerivIn2 >
  409. value_type prepare_dense_output( int k , const StateIn1 &x_start , const DerivIn1 &dxdt_start ,
  410. const StateIn2 & /* x_end */ , const DerivIn2 & /*dxdt_end */ , time_type dt )
  411. /* k is the order to which the result was approximated */
  412. {
  413. /* compute the coefficients of the interpolation polynomial
  414. * we parametrize the interval t .. t+dt by theta = -1 .. 1
  415. * we use 2k+3 values at the interval center theta=0 to obtain the interpolation coefficients
  416. * the values are x(t+dt/2) and the derivatives dx/dt , ... d^(2k+2) x / dt^(2k+2) at the midpoints
  417. * the derivatives are approximated via finite differences
  418. * all values are obtained from interpolation of the results from the increasing orders of the midpoint calls
  419. */
  420. // calculate finite difference approximations to derivatives at the midpoint
  421. for( int j = 0 ; j<=k ; j++ )
  422. {
  423. /* not working with boost units... */
  424. const value_type d = m_interval_sequence[j] / ( static_cast<value_type>(2) * dt );
  425. value_type f = 1.0; //factor 1/2 here because our interpolation interval has length 2 !!!
  426. for( int kappa = 0 ; kappa <= 2*j+1 ; ++kappa )
  427. {
  428. calculate_finite_difference( j , kappa , f , dxdt_start );
  429. f *= d;
  430. }
  431. if( j > 0 )
  432. extrapolate_dense_out( j , m_mp_states , m_coeff );
  433. }
  434. time_type d = dt/2;
  435. // extrapolate finite differences
  436. for( int kappa = 0 ; kappa<=2*k+1 ; kappa++ )
  437. {
  438. for( int j=1 ; j<=(k-kappa/2) ; ++j )
  439. extrapolate_dense_out( j , m_diffs[kappa] , m_coeff , kappa/2 );
  440. // extrapolation results are now stored in m_diffs[kappa][0]
  441. // divide kappa-th derivative by kappa because we need these terms for dense output interpolation
  442. m_algebra.for_each1( m_diffs[kappa][0].m_v , typename operations_type::template scale< time_type >( static_cast<time_type>(d) ) );
  443. d *= dt/(2*(kappa+2));
  444. }
  445. // dense output coefficients a_0 is stored in m_mp_states[0], a_i for i = 1...2k are stored in m_diffs[i-1][0]
  446. // the error is just the highest order coefficient of the interpolation polynomial
  447. // this is because we use only the midpoint theta=0 as support for the interpolation (remember that theta = -1 .. 1)
  448. value_type error = 0.0;
  449. if( m_control_interpolation )
  450. {
  451. boost::numeric::odeint::copy( m_diffs[2*k+1][0].m_v , m_err.m_v );
  452. error = m_error_checker.error( m_algebra , x_start , dxdt_start , m_err.m_v , dt );
  453. }
  454. return error;
  455. }
  456. template< class DerivIn >
  457. void calculate_finite_difference( size_t j , size_t kappa , value_type fac , const DerivIn &dxdt )
  458. {
  459. const int m = m_interval_sequence[j]/2-1;
  460. if( kappa == 0) // no calculation required for 0th derivative of f
  461. {
  462. m_algebra.for_each2( m_diffs[0][j].m_v , m_derivs[j][m].m_v ,
  463. typename operations_type::template scale_sum1< value_type >( fac ) );
  464. }
  465. else
  466. {
  467. // calculate the index of m_diffs for this kappa-j-combination
  468. const int j_diffs = j - kappa/2;
  469. m_algebra.for_each2( m_diffs[kappa][j_diffs].m_v , m_derivs[j][m+kappa].m_v ,
  470. typename operations_type::template scale_sum1< value_type >( fac ) );
  471. value_type sign = -1.0;
  472. int c = 1;
  473. //computes the j-th order finite difference for the kappa-th derivative of f at t+dt/2 using function evaluations stored in m_derivs
  474. for( int i = m+static_cast<int>(kappa)-2 ; i >= m-static_cast<int>(kappa) ; i -= 2 )
  475. {
  476. if( i >= 0 )
  477. {
  478. m_algebra.for_each3( m_diffs[kappa][j_diffs].m_v , m_diffs[kappa][j_diffs].m_v , m_derivs[j][i].m_v ,
  479. typename operations_type::template scale_sum2< value_type , value_type >( 1.0 ,
  480. sign * fac * boost::math::binomial_coefficient< value_type >( kappa , c ) ) );
  481. }
  482. else
  483. {
  484. m_algebra.for_each3( m_diffs[kappa][j_diffs].m_v , m_diffs[kappa][j_diffs].m_v , dxdt ,
  485. typename operations_type::template scale_sum2< value_type , value_type >( 1.0 , sign * fac ) );
  486. }
  487. sign *= -1;
  488. ++c;
  489. }
  490. }
  491. }
  492. template< class StateOut >
  493. void do_interpolation( time_type t , StateOut &out ) const
  494. {
  495. // interpolation polynomial is defined for theta = -1 ... 1
  496. // m_k_final is the number of order-iterations done for the last step - it governs the order of the interpolation polynomial
  497. const value_type theta = 2 * get_unit_value( (t - m_t_last) / (m_t - m_t_last) ) - 1;
  498. // we use only values at interval center, that is theta=0, for interpolation
  499. // our interpolation polynomial is thus of order 2k+2, hence we have 2k+3 terms
  500. boost::numeric::odeint::copy( m_mp_states[0].m_v , out );
  501. // add remaining terms: x += a_1 theta + a2 theta^2 + ... + a_{2k} theta^{2k}
  502. value_type theta_pow( theta );
  503. for( size_t i=0 ; i<=2*m_k_final+1 ; ++i )
  504. {
  505. m_algebra.for_each3( out , out , m_diffs[i][0].m_v ,
  506. typename operations_type::template scale_sum2< value_type >( static_cast<value_type>(1) , theta_pow ) );
  507. theta_pow *= theta;
  508. }
  509. }
  510. /* Resizer methods */
  511. template< class StateIn >
  512. bool resize_impl( const StateIn &x )
  513. {
  514. bool resized( false );
  515. resized |= adjust_size_by_resizeability( m_x1 , x , typename is_resizeable<state_type>::type() );
  516. resized |= adjust_size_by_resizeability( m_x2 , x , typename is_resizeable<state_type>::type() );
  517. resized |= adjust_size_by_resizeability( m_dxdt1 , x , typename is_resizeable<state_type>::type() );
  518. resized |= adjust_size_by_resizeability( m_dxdt2 , x , typename is_resizeable<state_type>::type() );
  519. resized |= adjust_size_by_resizeability( m_err , x , typename is_resizeable<state_type>::type() );
  520. for( size_t i = 0 ; i < m_k_max ; ++i )
  521. resized |= adjust_size_by_resizeability( m_table[i] , x , typename is_resizeable<state_type>::type() );
  522. for( size_t i = 0 ; i < m_k_max+1 ; ++i )
  523. resized |= adjust_size_by_resizeability( m_mp_states[i] , x , typename is_resizeable<state_type>::type() );
  524. for( size_t i = 0 ; i < m_k_max+1 ; ++i )
  525. for( size_t j = 0 ; j < m_derivs[i].size() ; ++j )
  526. resized |= adjust_size_by_resizeability( m_derivs[i][j] , x , typename is_resizeable<deriv_type>::type() );
  527. for( size_t i = 0 ; i < 2*m_k_max+2 ; ++i )
  528. for( size_t j = 0 ; j < m_diffs[i].size() ; ++j )
  529. resized |= adjust_size_by_resizeability( m_diffs[i][j] , x , typename is_resizeable<deriv_type>::type() );
  530. return resized;
  531. }
  532. state_type& get_current_state( void )
  533. {
  534. return m_current_state_x1 ? m_x1.m_v : m_x2.m_v ;
  535. }
  536. const state_type& get_current_state( void ) const
  537. {
  538. return m_current_state_x1 ? m_x1.m_v : m_x2.m_v ;
  539. }
  540. state_type& get_old_state( void )
  541. {
  542. return m_current_state_x1 ? m_x2.m_v : m_x1.m_v ;
  543. }
  544. const state_type& get_old_state( void ) const
  545. {
  546. return m_current_state_x1 ? m_x2.m_v : m_x1.m_v ;
  547. }
  548. deriv_type& get_current_deriv( void )
  549. {
  550. return m_current_state_x1 ? m_dxdt1.m_v : m_dxdt2.m_v ;
  551. }
  552. const deriv_type& get_current_deriv( void ) const
  553. {
  554. return m_current_state_x1 ? m_dxdt1.m_v : m_dxdt2.m_v ;
  555. }
  556. deriv_type& get_old_deriv( void )
  557. {
  558. return m_current_state_x1 ? m_dxdt2.m_v : m_dxdt1.m_v ;
  559. }
  560. const deriv_type& get_old_deriv( void ) const
  561. {
  562. return m_current_state_x1 ? m_dxdt2.m_v : m_dxdt1.m_v ;
  563. }
  564. void toggle_current_state( void )
  565. {
  566. m_current_state_x1 = ! m_current_state_x1;
  567. }
  568. default_error_checker< value_type, algebra_type , operations_type > m_error_checker;
  569. modified_midpoint_dense_out< state_type , value_type , deriv_type , time_type , algebra_type , operations_type , resizer_type > m_midpoint;
  570. bool m_control_interpolation;
  571. bool m_last_step_rejected;
  572. bool m_first;
  573. time_type m_t;
  574. time_type m_dt;
  575. time_type m_dt_last;
  576. time_type m_t_last;
  577. size_t m_current_k_opt;
  578. size_t m_k_final;
  579. algebra_type m_algebra;
  580. resizer_type m_resizer;
  581. wrapped_state_type m_x1 , m_x2;
  582. wrapped_deriv_type m_dxdt1 , m_dxdt2;
  583. wrapped_state_type m_err;
  584. bool m_current_state_x1;
  585. value_vector m_error; // errors of repeated midpoint steps and extrapolations
  586. int_vector m_interval_sequence; // stores the successive interval counts
  587. value_matrix m_coeff;
  588. int_vector m_cost; // costs for interval count
  589. value_vector m_facmin_table; // for precomputed facmin to save pow calls
  590. state_vector_type m_table; // sequence of states for extrapolation
  591. //for dense output:
  592. state_vector_type m_mp_states; // sequence of approximations of x at distance center
  593. deriv_table_type m_derivs; // table of function values
  594. deriv_table_type m_diffs; // table of function values
  595. //wrapped_state_type m_a1 , m_a2 , m_a3 , m_a4;
  596. value_type STEPFAC1 , STEPFAC2 , STEPFAC3 , STEPFAC4 , KFAC1 , KFAC2;
  597. };
  598. /********** DOXYGEN **********/
  599. /**
  600. * \class bulirsch_stoer_dense_out
  601. * \brief The Bulirsch-Stoer algorithm.
  602. *
  603. * The Bulirsch-Stoer is a controlled stepper that adjusts both step size
  604. * and order of the method. The algorithm uses the modified midpoint and
  605. * a polynomial extrapolation compute the solution. This class also provides
  606. * dense output facility.
  607. *
  608. * \tparam State The state type.
  609. * \tparam Value The value type.
  610. * \tparam Deriv The type representing the time derivative of the state.
  611. * \tparam Time The time representing the independent variable - the time.
  612. * \tparam Algebra The algebra type.
  613. * \tparam Operations The operations type.
  614. * \tparam Resizer The resizer policy type.
  615. */
  616. /**
  617. * \fn bulirsch_stoer_dense_out::bulirsch_stoer_dense_out( value_type eps_abs , value_type eps_rel , value_type factor_x , value_type factor_dxdt , bool control_interpolation )
  618. * \brief Constructs the bulirsch_stoer class, including initialization of
  619. * the error bounds.
  620. *
  621. * \param eps_abs Absolute tolerance level.
  622. * \param eps_rel Relative tolerance level.
  623. * \param factor_x Factor for the weight of the state.
  624. * \param factor_dxdt Factor for the weight of the derivative.
  625. * \param control_interpolation Set true to additionally control the error of
  626. * the interpolation.
  627. */
  628. /**
  629. * \fn bulirsch_stoer_dense_out::try_step( System system , const StateIn &in , const DerivIn &dxdt , time_type &t , StateOut &out , DerivOut &dxdt_new , time_type &dt )
  630. * \brief Tries to perform one step.
  631. *
  632. * This method tries to do one step with step size dt. If the error estimate
  633. * is to large, the step is rejected and the method returns fail and the
  634. * step size dt is reduced. If the error estimate is acceptably small, the
  635. * step is performed, success is returned and dt might be increased to make
  636. * the steps as large as possible. This method also updates t if a step is
  637. * performed. Also, the internal order of the stepper is adjusted if required.
  638. *
  639. * \param system The system function to solve, hence the r.h.s. of the ODE.
  640. * It must fulfill the Simple System concept.
  641. * \param in The state of the ODE which should be solved.
  642. * \param dxdt The derivative of state.
  643. * \param t The value of the time. Updated if the step is successful.
  644. * \param out Used to store the result of the step.
  645. * \param dt The step size. Updated.
  646. * \return success if the step was accepted, fail otherwise.
  647. */
  648. /**
  649. * \fn bulirsch_stoer_dense_out::initialize( const StateType &x0 , const time_type &t0 , const time_type &dt0 )
  650. * \brief Initializes the dense output stepper.
  651. *
  652. * \param x0 The initial state.
  653. * \param t0 The initial time.
  654. * \param dt0 The initial time step.
  655. */
  656. /**
  657. * \fn bulirsch_stoer_dense_out::do_step( System system )
  658. * \brief Does one time step. This is the main method that should be used to
  659. * integrate an ODE with this stepper.
  660. * \note initialize has to be called before using this method to set the
  661. * initial conditions x,t and the stepsize.
  662. * \param system The system function to solve, hence the r.h.s. of the
  663. * ordinary differential equation. It must fulfill the Simple System concept.
  664. * \return Pair with start and end time of the integration step.
  665. */
  666. /**
  667. * \fn bulirsch_stoer_dense_out::calc_state( time_type t , StateOut &x ) const
  668. * \brief Calculates the solution at an intermediate point within the last step
  669. * \param t The time at which the solution should be calculated, has to be
  670. * in the current time interval.
  671. * \param x The output variable where the result is written into.
  672. */
  673. /**
  674. * \fn bulirsch_stoer_dense_out::current_state( void ) const
  675. * \brief Returns the current state of the solution.
  676. * \return The current state of the solution x(t).
  677. */
  678. /**
  679. * \fn bulirsch_stoer_dense_out::current_time( void ) const
  680. * \brief Returns the current time of the solution.
  681. * \return The current time of the solution t.
  682. */
  683. /**
  684. * \fn bulirsch_stoer_dense_out::previous_state( void ) const
  685. * \brief Returns the last state of the solution.
  686. * \return The last state of the solution x(t-dt).
  687. */
  688. /**
  689. * \fn bulirsch_stoer_dense_out::previous_time( void ) const
  690. * \brief Returns the last time of the solution.
  691. * \return The last time of the solution t-dt.
  692. */
  693. /**
  694. * \fn bulirsch_stoer_dense_out::current_time_step( void ) const
  695. * \brief Returns the current step size.
  696. * \return The current step size.
  697. */
  698. /**
  699. * \fn bulirsch_stoer_dense_out::adjust_size( const StateIn &x )
  700. * \brief Adjust the size of all temporaries in the stepper manually.
  701. * \param x A state from which the size of the temporaries to be resized is deduced.
  702. */
  703. }
  704. }
  705. }
  706. #endif // BOOST_NUMERIC_ODEINT_STEPPER_BULIRSCH_STOER_HPP_INCLUDED