bulirsch_stoer.hpp 26 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642
  1. /*
  2. [auto_generated]
  3. boost/numeric/odeint/stepper/bulirsch_stoer.hpp
  4. [begin_description]
  5. Implementation of the Burlish-Stoer method. As described in
  6. Ernst Hairer, Syvert Paul Norsett, Gerhard Wanner
  7. Solving Ordinary Differential Equations I. Nonstiff Problems.
  8. Springer Series in Comput. Mathematics, Vol. 8, Springer-Verlag 1987, Second revised edition 1993.
  9. [end_description]
  10. Copyright 2011-2013 Mario Mulansky
  11. Copyright 2011-2013 Karsten Ahnert
  12. Copyright 2012 Christoph Koke
  13. Distributed under the Boost Software License, Version 1.0.
  14. (See accompanying file LICENSE_1_0.txt or
  15. copy at http://www.boost.org/LICENSE_1_0.txt)
  16. */
  17. #ifndef BOOST_NUMERIC_ODEINT_STEPPER_BULIRSCH_STOER_HPP_INCLUDED
  18. #define BOOST_NUMERIC_ODEINT_STEPPER_BULIRSCH_STOER_HPP_INCLUDED
  19. #include <iostream>
  20. #include <algorithm>
  21. #include <boost/config.hpp> // for min/max guidelines
  22. #include <boost/numeric/odeint/util/bind.hpp>
  23. #include <boost/numeric/odeint/util/unwrap_reference.hpp>
  24. #include <boost/numeric/odeint/stepper/controlled_runge_kutta.hpp>
  25. #include <boost/numeric/odeint/stepper/modified_midpoint.hpp>
  26. #include <boost/numeric/odeint/stepper/controlled_step_result.hpp>
  27. #include <boost/numeric/odeint/algebra/range_algebra.hpp>
  28. #include <boost/numeric/odeint/algebra/default_operations.hpp>
  29. #include <boost/numeric/odeint/algebra/algebra_dispatcher.hpp>
  30. #include <boost/numeric/odeint/algebra/operations_dispatcher.hpp>
  31. #include <boost/numeric/odeint/util/state_wrapper.hpp>
  32. #include <boost/numeric/odeint/util/is_resizeable.hpp>
  33. #include <boost/numeric/odeint/util/resizer.hpp>
  34. #include <boost/numeric/odeint/util/unit_helper.hpp>
  35. #include <boost/numeric/odeint/util/detail/less_with_sign.hpp>
  36. namespace boost {
  37. namespace numeric {
  38. namespace odeint {
  39. template<
  40. class State ,
  41. class Value = double ,
  42. class Deriv = State ,
  43. class Time = Value ,
  44. class Algebra = typename algebra_dispatcher< State >::algebra_type ,
  45. class Operations = typename operations_dispatcher< State >::operations_type ,
  46. class Resizer = initially_resizer
  47. >
  48. class bulirsch_stoer {
  49. public:
  50. typedef State state_type;
  51. typedef Value value_type;
  52. typedef Deriv deriv_type;
  53. typedef Time time_type;
  54. typedef Algebra algebra_type;
  55. typedef Operations operations_type;
  56. typedef Resizer resizer_type;
  57. #ifndef DOXYGEN_SKIP
  58. typedef state_wrapper< state_type > wrapped_state_type;
  59. typedef state_wrapper< deriv_type > wrapped_deriv_type;
  60. typedef controlled_stepper_tag stepper_category;
  61. typedef bulirsch_stoer< State , Value , Deriv , Time , Algebra , Operations , Resizer > controlled_error_bs_type;
  62. typedef typename inverse_time< time_type >::type inv_time_type;
  63. typedef std::vector< value_type > value_vector;
  64. typedef std::vector< time_type > time_vector;
  65. typedef std::vector< inv_time_type > inv_time_vector; //should be 1/time_type for boost.units
  66. typedef std::vector< value_vector > value_matrix;
  67. typedef std::vector< size_t > int_vector;
  68. typedef std::vector< wrapped_state_type > state_table_type;
  69. #endif //DOXYGEN_SKIP
  70. const static size_t m_k_max = 8;
  71. bulirsch_stoer(
  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. : m_error_checker( eps_abs , eps_rel , factor_x, factor_dxdt ) , m_midpoint() ,
  76. m_last_step_rejected( false ) , m_first( true ) ,
  77. m_max_dt(max_dt) ,
  78. m_interval_sequence( m_k_max+1 ) ,
  79. m_coeff( m_k_max+1 ) ,
  80. m_cost( m_k_max+1 ) ,
  81. m_facmin_table( m_k_max+1 ) ,
  82. m_table( m_k_max ) ,
  83. STEPFAC1( 0.65 ) , STEPFAC2( 0.94 ) , STEPFAC3( 0.02 ) , STEPFAC4( 4.0 ) , KFAC1( 0.8 ) , KFAC2( 0.9 )
  84. {
  85. BOOST_USING_STD_MIN();
  86. BOOST_USING_STD_MAX();
  87. /* initialize sequence of stage numbers and work */
  88. for( unsigned short i = 0; i < m_k_max+1; i++ )
  89. {
  90. m_interval_sequence[i] = 2 * (i+1);
  91. if( i == 0 )
  92. m_cost[i] = m_interval_sequence[i];
  93. else
  94. m_cost[i] = m_cost[i-1] + m_interval_sequence[i];
  95. m_coeff[i].resize(i);
  96. m_facmin_table[i] = pow BOOST_PREVENT_MACRO_SUBSTITUTION( STEPFAC3 , static_cast< value_type >(1) / static_cast< value_type >( 2*i+1 ) );
  97. for( size_t k = 0 ; k < i ; ++k )
  98. {
  99. const value_type r = static_cast< value_type >( m_interval_sequence[i] ) / static_cast< value_type >( m_interval_sequence[k] );
  100. m_coeff[i][k] = 1.0 / ( r*r - static_cast< value_type >( 1.0 ) ); // coefficients for extrapolation
  101. }
  102. }
  103. reset();
  104. }
  105. /*
  106. * Version 1 : try_step( sys , x , t , dt )
  107. *
  108. * The overloads are needed to solve the forwarding problem
  109. */
  110. template< class System , class StateInOut >
  111. controlled_step_result try_step( System system , StateInOut &x , time_type &t , time_type &dt )
  112. {
  113. return try_step_v1( system , x , t, dt );
  114. }
  115. /**
  116. * \brief Second version to solve the forwarding problem, can be used with Boost.Range as StateInOut.
  117. */
  118. template< class System , class StateInOut >
  119. controlled_step_result try_step( System system , const StateInOut &x , time_type &t , time_type &dt )
  120. {
  121. return try_step_v1( system , x , t, dt );
  122. }
  123. /*
  124. * Version 2 : try_step( sys , x , dxdt , t , dt )
  125. *
  126. * this version does not solve the forwarding problem, boost.range can not be used
  127. */
  128. template< class System , class StateInOut , class DerivIn >
  129. controlled_step_result try_step( System system , StateInOut &x , const DerivIn &dxdt , time_type &t , time_type &dt )
  130. {
  131. m_xnew_resizer.adjust_size( x , detail::bind( &controlled_error_bs_type::template resize_m_xnew< StateInOut > , detail::ref( *this ) , detail::_1 ) );
  132. controlled_step_result res = try_step( system , x , dxdt , t , m_xnew.m_v , dt );
  133. if( res == success )
  134. {
  135. boost::numeric::odeint::copy( m_xnew.m_v , x );
  136. }
  137. return res;
  138. }
  139. /*
  140. * Version 3 : try_step( sys , in , t , out , dt )
  141. *
  142. * this version does not solve the forwarding problem, boost.range can not be used
  143. */
  144. template< class System , class StateIn , class StateOut >
  145. typename boost::disable_if< boost::is_same< StateIn , time_type > , controlled_step_result >::type
  146. try_step( System system , const StateIn &in , time_type &t , StateOut &out , time_type &dt )
  147. {
  148. typename odeint::unwrap_reference< System >::type &sys = system;
  149. m_dxdt_resizer.adjust_size( in , detail::bind( &controlled_error_bs_type::template resize_m_dxdt< StateIn > , detail::ref( *this ) , detail::_1 ) );
  150. sys( in , m_dxdt.m_v , t );
  151. return try_step( system , in , m_dxdt.m_v , t , out , dt );
  152. }
  153. /*
  154. * Full version : try_step( sys , in , dxdt_in , t , out , dt )
  155. *
  156. * contains the actual implementation
  157. */
  158. template< class System , class StateIn , class DerivIn , class StateOut >
  159. controlled_step_result try_step( System system , const StateIn &in , const DerivIn &dxdt , time_type &t , StateOut &out , time_type &dt )
  160. {
  161. if( m_max_dt != static_cast<time_type>(0) && detail::less_with_sign(m_max_dt, dt, dt) )
  162. {
  163. // given step size is bigger then max_dt
  164. // set limit and return fail
  165. dt = m_max_dt;
  166. return fail;
  167. }
  168. BOOST_USING_STD_MIN();
  169. BOOST_USING_STD_MAX();
  170. static const value_type val1( 1.0 );
  171. if( m_resizer.adjust_size( in , detail::bind( &controlled_error_bs_type::template resize_impl< StateIn > , detail::ref( *this ) , detail::_1 ) ) )
  172. {
  173. reset(); // system resized -> reset
  174. }
  175. if( dt != m_dt_last )
  176. {
  177. reset(); // step size changed from outside -> reset
  178. }
  179. bool reject( true );
  180. time_vector h_opt( m_k_max+1 );
  181. inv_time_vector work( m_k_max+1 );
  182. time_type new_h = dt;
  183. /* m_current_k_opt is the estimated current optimal stage number */
  184. for( size_t k = 0 ; k <= m_current_k_opt+1 ; k++ )
  185. {
  186. /* the stage counts are stored in m_interval_sequence */
  187. m_midpoint.set_steps( m_interval_sequence[k] );
  188. if( k == 0 )
  189. {
  190. m_midpoint.do_step( system , in , dxdt , t , out , dt );
  191. /* the first step, nothing more to do */
  192. }
  193. else
  194. {
  195. m_midpoint.do_step( system , in , dxdt , t , m_table[k-1].m_v , dt );
  196. extrapolate( k , m_table , m_coeff , out );
  197. // get error estimate
  198. m_algebra.for_each3( m_err.m_v , out , m_table[0].m_v ,
  199. typename operations_type::template scale_sum2< value_type , value_type >( val1 , -val1 ) );
  200. const value_type error = m_error_checker.error( m_algebra , in , dxdt , m_err.m_v , dt );
  201. h_opt[k] = calc_h_opt( dt , error , k );
  202. work[k] = static_cast<value_type>( m_cost[k] ) / h_opt[k];
  203. if( (k == m_current_k_opt-1) || m_first )
  204. { // convergence before k_opt ?
  205. if( error < 1.0 )
  206. {
  207. //convergence
  208. reject = false;
  209. if( (work[k] < KFAC2*work[k-1]) || (m_current_k_opt <= 2) )
  210. {
  211. // leave order as is (except we were in first round)
  212. 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 ) );
  213. new_h = h_opt[k];
  214. new_h *= static_cast<value_type>( m_cost[k+1] ) / static_cast<value_type>( m_cost[k] );
  215. } else {
  216. 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) ) );
  217. new_h = h_opt[k];
  218. }
  219. break;
  220. }
  221. else if( should_reject( error , k ) && !m_first )
  222. {
  223. reject = true;
  224. new_h = h_opt[k];
  225. break;
  226. }
  227. }
  228. if( k == m_current_k_opt )
  229. { // convergence at k_opt ?
  230. if( error < 1.0 )
  231. {
  232. //convergence
  233. reject = false;
  234. if( (work[k-1] < KFAC2*work[k]) )
  235. {
  236. m_current_k_opt = max BOOST_PREVENT_MACRO_SUBSTITUTION( 2 , static_cast<int>(m_current_k_opt)-1 );
  237. new_h = h_opt[m_current_k_opt];
  238. }
  239. else if( (work[k] < KFAC2*work[k-1]) && !m_last_step_rejected )
  240. {
  241. m_current_k_opt = min BOOST_PREVENT_MACRO_SUBSTITUTION( static_cast<int>(m_k_max-1) , static_cast<int>(m_current_k_opt)+1 );
  242. new_h = h_opt[k];
  243. new_h *= static_cast<value_type>(m_cost[m_current_k_opt])/static_cast<value_type>(m_cost[k]);
  244. } else
  245. new_h = h_opt[m_current_k_opt];
  246. break;
  247. }
  248. else if( should_reject( error , k ) )
  249. {
  250. reject = true;
  251. new_h = h_opt[m_current_k_opt];
  252. break;
  253. }
  254. }
  255. if( k == m_current_k_opt+1 )
  256. { // convergence at k_opt+1 ?
  257. if( error < 1.0 )
  258. { //convergence
  259. reject = false;
  260. if( work[k-2] < KFAC2*work[k-1] )
  261. m_current_k_opt = max BOOST_PREVENT_MACRO_SUBSTITUTION( 2 , static_cast<int>(m_current_k_opt)-1 );
  262. if( (work[k] < KFAC2*work[m_current_k_opt]) && !m_last_step_rejected )
  263. m_current_k_opt = min BOOST_PREVENT_MACRO_SUBSTITUTION( static_cast<int>(m_k_max)-1 , static_cast<int>(k) );
  264. new_h = h_opt[m_current_k_opt];
  265. } else
  266. {
  267. reject = true;
  268. new_h = h_opt[m_current_k_opt];
  269. }
  270. break;
  271. }
  272. }
  273. }
  274. if( !reject )
  275. {
  276. t += dt;
  277. }
  278. if( !m_last_step_rejected || boost::numeric::odeint::detail::less_with_sign(new_h, dt, dt) )
  279. {
  280. // limit step size
  281. if( m_max_dt != static_cast<time_type>(0) )
  282. {
  283. new_h = detail::min_abs(m_max_dt, new_h);
  284. }
  285. m_dt_last = new_h;
  286. dt = new_h;
  287. }
  288. m_last_step_rejected = reject;
  289. m_first = false;
  290. if( reject )
  291. return fail;
  292. else
  293. return success;
  294. }
  295. /** \brief Resets the internal state of the stepper */
  296. void reset()
  297. {
  298. m_first = true;
  299. m_last_step_rejected = false;
  300. // crude estimate of optimal order
  301. m_current_k_opt = 4;
  302. /* no calculation because log10 might not exist for value_type!
  303. const value_type logfact( -log10( max BOOST_PREVENT_MACRO_SUBSTITUTION( eps_rel , static_cast< value_type >(1.0E-12) ) ) * 0.6 + 0.5 );
  304. m_current_k_opt = max BOOST_PREVENT_MACRO_SUBSTITUTION( static_cast<value_type>( 1 ) , min BOOST_PREVENT_MACRO_SUBSTITUTION( static_cast<value_type>( m_k_max-1 ) , logfact ));
  305. */
  306. }
  307. /* Resizer methods */
  308. template< class StateIn >
  309. void adjust_size( const StateIn &x )
  310. {
  311. resize_m_dxdt( x );
  312. resize_m_xnew( x );
  313. resize_impl( x );
  314. m_midpoint.adjust_size( x );
  315. }
  316. private:
  317. template< class StateIn >
  318. bool resize_m_dxdt( const StateIn &x )
  319. {
  320. return adjust_size_by_resizeability( m_dxdt , x , typename is_resizeable<deriv_type>::type() );
  321. }
  322. template< class StateIn >
  323. bool resize_m_xnew( const StateIn &x )
  324. {
  325. return adjust_size_by_resizeability( m_xnew , x , typename is_resizeable<state_type>::type() );
  326. }
  327. template< class StateIn >
  328. bool resize_impl( const StateIn &x )
  329. {
  330. bool resized( false );
  331. for( size_t i = 0 ; i < m_k_max ; ++i )
  332. resized |= adjust_size_by_resizeability( m_table[i] , x , typename is_resizeable<state_type>::type() );
  333. resized |= adjust_size_by_resizeability( m_err , x , typename is_resizeable<state_type>::type() );
  334. return resized;
  335. }
  336. template< class System , class StateInOut >
  337. controlled_step_result try_step_v1( System system , StateInOut &x , time_type &t , time_type &dt )
  338. {
  339. typename odeint::unwrap_reference< System >::type &sys = system;
  340. m_dxdt_resizer.adjust_size( x , detail::bind( &controlled_error_bs_type::template resize_m_dxdt< StateInOut > , detail::ref( *this ) , detail::_1 ) );
  341. sys( x , m_dxdt.m_v ,t );
  342. return try_step( system , x , m_dxdt.m_v , t , dt );
  343. }
  344. template< class StateInOut >
  345. void extrapolate( size_t k , state_table_type &table , const value_matrix &coeff , StateInOut &xest )
  346. /* polynomial extrapolation, see http://www.nr.com/webnotes/nr3web21.pdf
  347. uses the obtained intermediate results to extrapolate to dt->0
  348. */
  349. {
  350. static const value_type val1 = static_cast< value_type >( 1.0 );
  351. for( int j=k-1 ; j>0 ; --j )
  352. {
  353. m_algebra.for_each3( table[j-1].m_v , table[j].m_v , table[j-1].m_v ,
  354. typename operations_type::template scale_sum2< value_type , value_type >( val1 + coeff[k][j] , -coeff[k][j] ) );
  355. }
  356. m_algebra.for_each3( xest , table[0].m_v , xest ,
  357. typename operations_type::template scale_sum2< value_type , value_type >( val1 + coeff[k][0] , -coeff[k][0]) );
  358. }
  359. time_type calc_h_opt( time_type h , value_type error , size_t k ) const
  360. /* calculates the optimal step size for a given error and stage number */
  361. {
  362. BOOST_USING_STD_MIN();
  363. BOOST_USING_STD_MAX();
  364. using std::pow;
  365. value_type expo( 1.0/(2*k+1) );
  366. value_type facmin = m_facmin_table[k];
  367. value_type fac;
  368. if (error == 0.0)
  369. fac=1.0/facmin;
  370. else
  371. {
  372. fac = STEPFAC2 / pow BOOST_PREVENT_MACRO_SUBSTITUTION( error / STEPFAC1 , expo );
  373. fac = max BOOST_PREVENT_MACRO_SUBSTITUTION( static_cast<value_type>(facmin/STEPFAC4) , min BOOST_PREVENT_MACRO_SUBSTITUTION( static_cast<value_type>(1.0/facmin) , fac ) );
  374. }
  375. return h*fac;
  376. }
  377. controlled_step_result set_k_opt( size_t k , const inv_time_vector &work , const time_vector &h_opt , time_type &dt )
  378. /* calculates the optimal stage number */
  379. {
  380. if( k == 1 )
  381. {
  382. m_current_k_opt = 2;
  383. return success;
  384. }
  385. if( (work[k-1] < KFAC1*work[k]) || (k == m_k_max) )
  386. { // order decrease
  387. m_current_k_opt = k-1;
  388. dt = h_opt[ m_current_k_opt ];
  389. return success;
  390. }
  391. else if( (work[k] < KFAC2*work[k-1]) || m_last_step_rejected || (k == m_k_max-1) )
  392. { // same order - also do this if last step got rejected
  393. m_current_k_opt = k;
  394. dt = h_opt[ m_current_k_opt ];
  395. return success;
  396. }
  397. else
  398. { // order increase - only if last step was not rejected
  399. m_current_k_opt = k+1;
  400. dt = h_opt[ m_current_k_opt-1 ] * m_cost[ m_current_k_opt ] / m_cost[ m_current_k_opt-1 ] ;
  401. return success;
  402. }
  403. }
  404. bool in_convergence_window( size_t k ) const
  405. {
  406. if( (k == m_current_k_opt-1) && !m_last_step_rejected )
  407. return true; // decrease stepsize only if last step was not rejected
  408. return ( (k == m_current_k_opt) || (k == m_current_k_opt+1) );
  409. }
  410. bool should_reject( value_type error , size_t k ) const
  411. {
  412. if( k == m_current_k_opt-1 )
  413. {
  414. const value_type d = m_interval_sequence[m_current_k_opt] * m_interval_sequence[m_current_k_opt+1] /
  415. (m_interval_sequence[0]*m_interval_sequence[0]);
  416. //step will fail, criterion 17.3.17 in NR
  417. return ( error > d*d );
  418. }
  419. else if( k == m_current_k_opt )
  420. {
  421. const value_type d = m_interval_sequence[m_current_k_opt] / m_interval_sequence[0];
  422. return ( error > d*d );
  423. } else
  424. return error > 1.0;
  425. }
  426. default_error_checker< value_type, algebra_type , operations_type > m_error_checker;
  427. modified_midpoint< state_type , value_type , deriv_type , time_type , algebra_type , operations_type , resizer_type > m_midpoint;
  428. bool m_last_step_rejected;
  429. bool m_first;
  430. time_type m_dt_last;
  431. time_type m_t_last;
  432. time_type m_max_dt;
  433. size_t m_current_k_opt;
  434. algebra_type m_algebra;
  435. resizer_type m_dxdt_resizer;
  436. resizer_type m_xnew_resizer;
  437. resizer_type m_resizer;
  438. wrapped_state_type m_xnew;
  439. wrapped_state_type m_err;
  440. wrapped_deriv_type m_dxdt;
  441. int_vector m_interval_sequence; // stores the successive interval counts
  442. value_matrix m_coeff;
  443. int_vector m_cost; // costs for interval count
  444. value_vector m_facmin_table; // for precomputed facmin to save pow calls
  445. state_table_type m_table; // sequence of states for extrapolation
  446. value_type STEPFAC1 , STEPFAC2 , STEPFAC3 , STEPFAC4 , KFAC1 , KFAC2;
  447. };
  448. /******** DOXYGEN ********/
  449. /**
  450. * \class bulirsch_stoer
  451. * \brief The Bulirsch-Stoer algorithm.
  452. *
  453. * The Bulirsch-Stoer is a controlled stepper that adjusts both step size
  454. * and order of the method. The algorithm uses the modified midpoint and
  455. * a polynomial extrapolation compute the solution.
  456. *
  457. * \tparam State The state type.
  458. * \tparam Value The value type.
  459. * \tparam Deriv The type representing the time derivative of the state.
  460. * \tparam Time The time representing the independent variable - the time.
  461. * \tparam Algebra The algebra type.
  462. * \tparam Operations The operations type.
  463. * \tparam Resizer The resizer policy type.
  464. */
  465. /**
  466. * \fn bulirsch_stoer::bulirsch_stoer( value_type eps_abs , value_type eps_rel , value_type factor_x , value_type factor_dxdt )
  467. * \brief Constructs the bulirsch_stoer class, including initialization of
  468. * the error bounds.
  469. *
  470. * \param eps_abs Absolute tolerance level.
  471. * \param eps_rel Relative tolerance level.
  472. * \param factor_x Factor for the weight of the state.
  473. * \param factor_dxdt Factor for the weight of the derivative.
  474. */
  475. /**
  476. * \fn bulirsch_stoer::try_step( System system , StateInOut &x , time_type &t , time_type &dt )
  477. * \brief Tries to perform one step.
  478. *
  479. * This method tries to do one step with step size dt. If the error estimate
  480. * is to large, the step is rejected and the method returns fail and the
  481. * step size dt is reduced. If the error estimate is acceptably small, the
  482. * step is performed, success is returned and dt might be increased to make
  483. * the steps as large as possible. This method also updates t if a step is
  484. * performed. Also, the internal order of the stepper is adjusted if required.
  485. *
  486. * \param system The system function to solve, hence the r.h.s. of the ODE.
  487. * It must fulfill the Simple System concept.
  488. * \param x The state of the ODE which should be solved. Overwritten if
  489. * the step is successful.
  490. * \param t The value of the time. Updated if the step is successful.
  491. * \param dt The step size. Updated.
  492. * \return success if the step was accepted, fail otherwise.
  493. */
  494. /**
  495. * \fn bulirsch_stoer::try_step( System system , StateInOut &x , const DerivIn &dxdt , time_type &t , time_type &dt )
  496. * \brief Tries to perform one step.
  497. *
  498. * This method tries to do one step with step size dt. If the error estimate
  499. * is to large, the step is rejected and the method returns fail and the
  500. * step size dt is reduced. If the error estimate is acceptably small, the
  501. * step is performed, success is returned and dt might be increased to make
  502. * the steps as large as possible. This method also updates t if a step is
  503. * performed. Also, the internal order of the stepper is adjusted if required.
  504. *
  505. * \param system The system function to solve, hence the r.h.s. of the ODE.
  506. * It must fulfill the Simple System concept.
  507. * \param x The state of the ODE which should be solved. Overwritten if
  508. * the step is successful.
  509. * \param dxdt The derivative of state.
  510. * \param t The value of the time. Updated if the step is successful.
  511. * \param dt The step size. Updated.
  512. * \return success if the step was accepted, fail otherwise.
  513. */
  514. /**
  515. * \fn bulirsch_stoer::try_step( System system , const StateIn &in , time_type &t , StateOut &out , time_type &dt )
  516. * \brief Tries to perform one step.
  517. *
  518. * \note This method is disabled if state_type=time_type to avoid ambiguity.
  519. *
  520. * This method tries to do one step with step size dt. If the error estimate
  521. * is to large, the step is rejected and the method returns fail and the
  522. * step size dt is reduced. If the error estimate is acceptably small, the
  523. * step is performed, success is returned and dt might be increased to make
  524. * the steps as large as possible. This method also updates t if a step is
  525. * performed. Also, the internal order of the stepper is adjusted if required.
  526. *
  527. * \param system The system function to solve, hence the r.h.s. of the ODE.
  528. * It must fulfill the Simple System concept.
  529. * \param in The state of the ODE which should be solved.
  530. * \param t The value of the time. Updated if the step is successful.
  531. * \param out Used to store the result of the step.
  532. * \param dt The step size. Updated.
  533. * \return success if the step was accepted, fail otherwise.
  534. */
  535. /**
  536. * \fn bulirsch_stoer::try_step( System system , const StateIn &in , const DerivIn &dxdt , time_type &t , StateOut &out , time_type &dt )
  537. * \brief Tries to perform one step.
  538. *
  539. * This method tries to do one step with step size dt. If the error estimate
  540. * is to large, the step is rejected and the method returns fail and the
  541. * step size dt is reduced. If the error estimate is acceptably small, the
  542. * step is performed, success is returned and dt might be increased to make
  543. * the steps as large as possible. This method also updates t if a step is
  544. * performed. Also, the internal order of the stepper is adjusted if required.
  545. *
  546. * \param system The system function to solve, hence the r.h.s. of the ODE.
  547. * It must fulfill the Simple System concept.
  548. * \param in The state of the ODE which should be solved.
  549. * \param dxdt The derivative of state.
  550. * \param t The value of the time. Updated if the step is successful.
  551. * \param out Used to store the result of the step.
  552. * \param dt The step size. Updated.
  553. * \return success if the step was accepted, fail otherwise.
  554. */
  555. /**
  556. * \fn bulirsch_stoer::adjust_size( const StateIn &x )
  557. * \brief Adjust the size of all temporaries in the stepper manually.
  558. * \param x A state from which the size of the temporaries to be resized is deduced.
  559. */
  560. }
  561. }
  562. }
  563. #endif // BOOST_NUMERIC_ODEINT_STEPPER_BULIRSCH_STOER_HPP_INCLUDED