bulirsch_stoer_dense_out.hpp 33 KB

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