controlled_runge_kutta.hpp 37 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018
  1. /* [auto_generated]
  2. boost/numeric/odeint/stepper/controlled_runge_kutta.hpp
  3. [begin_description]
  4. The default controlled stepper which can be used with all explicit Runge-Kutta error steppers.
  5. [end_description]
  6. Copyright 2010-2013 Karsten Ahnert
  7. Copyright 2010-2015 Mario Mulansky
  8. Copyright 2012 Christoph Koke
  9. Distributed under the Boost Software License, Version 1.0.
  10. (See accompanying file LICENSE_1_0.txt or
  11. copy at http://www.boost.org/LICENSE_1_0.txt)
  12. */
  13. #ifndef BOOST_NUMERIC_ODEINT_STEPPER_CONTROLLED_RUNGE_KUTTA_HPP_INCLUDED
  14. #define BOOST_NUMERIC_ODEINT_STEPPER_CONTROLLED_RUNGE_KUTTA_HPP_INCLUDED
  15. #include <cmath>
  16. #include <boost/config.hpp>
  17. #include <boost/utility/enable_if.hpp>
  18. #include <boost/type_traits/is_same.hpp>
  19. #include <boost/numeric/odeint/util/bind.hpp>
  20. #include <boost/numeric/odeint/util/unwrap_reference.hpp>
  21. #include <boost/numeric/odeint/util/copy.hpp>
  22. #include <boost/numeric/odeint/util/state_wrapper.hpp>
  23. #include <boost/numeric/odeint/util/is_resizeable.hpp>
  24. #include <boost/numeric/odeint/util/resizer.hpp>
  25. #include <boost/numeric/odeint/util/detail/less_with_sign.hpp>
  26. #include <boost/numeric/odeint/algebra/range_algebra.hpp>
  27. #include <boost/numeric/odeint/algebra/default_operations.hpp>
  28. #include <boost/numeric/odeint/algebra/algebra_dispatcher.hpp>
  29. #include <boost/numeric/odeint/stepper/controlled_step_result.hpp>
  30. #include <boost/numeric/odeint/stepper/stepper_categories.hpp>
  31. namespace boost {
  32. namespace numeric {
  33. namespace odeint {
  34. template
  35. <
  36. class Value ,
  37. class Algebra ,
  38. class Operations
  39. >
  40. class default_error_checker
  41. {
  42. public:
  43. typedef Value value_type;
  44. typedef Algebra algebra_type;
  45. typedef Operations operations_type;
  46. default_error_checker(
  47. value_type eps_abs = static_cast< value_type >( 1.0e-6 ) ,
  48. value_type eps_rel = static_cast< value_type >( 1.0e-6 ) ,
  49. value_type a_x = static_cast< value_type >( 1 ) ,
  50. value_type a_dxdt = static_cast< value_type >( 1 ))
  51. : m_eps_abs( eps_abs ) , m_eps_rel( eps_rel ) , m_a_x( a_x ) , m_a_dxdt( a_dxdt )
  52. { }
  53. template< class State , class Deriv , class Err, class Time >
  54. value_type error( const State &x_old , const Deriv &dxdt_old , Err &x_err , Time dt ) const
  55. {
  56. return error( algebra_type() , x_old , dxdt_old , x_err , dt );
  57. }
  58. template< class State , class Deriv , class Err, class Time >
  59. value_type error( algebra_type &algebra , const State &x_old , const Deriv &dxdt_old , Err &x_err , Time dt ) const
  60. {
  61. using std::abs;
  62. // this overwrites x_err !
  63. algebra.for_each3( x_err , x_old , dxdt_old ,
  64. typename operations_type::template rel_error< value_type >( m_eps_abs , m_eps_rel , m_a_x , m_a_dxdt * abs(get_unit_value( dt )) ) );
  65. // value_type res = algebra.reduce( x_err ,
  66. // typename operations_type::template maximum< value_type >() , static_cast< value_type >( 0 ) );
  67. return algebra.norm_inf( x_err );
  68. }
  69. private:
  70. value_type m_eps_abs;
  71. value_type m_eps_rel;
  72. value_type m_a_x;
  73. value_type m_a_dxdt;
  74. };
  75. template< typename Value, typename Time >
  76. class default_step_adjuster
  77. {
  78. public:
  79. typedef Time time_type;
  80. typedef Value value_type;
  81. default_step_adjuster(const time_type max_dt=static_cast<time_type>(0))
  82. : m_max_dt(max_dt)
  83. {}
  84. time_type decrease_step(time_type dt, const value_type error, const int error_order) const
  85. {
  86. // returns the decreased time step
  87. BOOST_USING_STD_MIN();
  88. BOOST_USING_STD_MAX();
  89. using std::pow;
  90. dt *= max
  91. BOOST_PREVENT_MACRO_SUBSTITUTION(
  92. static_cast<value_type>( static_cast<value_type>(9) / static_cast<value_type>(10) *
  93. pow(error, static_cast<value_type>(-1) / (error_order - 1))),
  94. static_cast<value_type>( static_cast<value_type>(1) / static_cast<value_type> (5)));
  95. if(m_max_dt != static_cast<time_type >(0))
  96. // limit to maximal stepsize even when decreasing
  97. dt = detail::min_abs(dt, m_max_dt);
  98. return dt;
  99. }
  100. time_type increase_step(time_type dt, value_type error, const int stepper_order) const
  101. {
  102. // returns the increased time step
  103. BOOST_USING_STD_MIN();
  104. BOOST_USING_STD_MAX();
  105. using std::pow;
  106. // adjust the size if dt is smaller than max_dt (providede max_dt is not zero)
  107. if(error < 0.5)
  108. {
  109. // error should be > 0
  110. error = max BOOST_PREVENT_MACRO_SUBSTITUTION (
  111. static_cast<value_type>( pow( static_cast<value_type>(5.0) , -static_cast<value_type>(stepper_order) ) ) ,
  112. error);
  113. // time_type dt_old = dt; unused variable warning
  114. //error too small - increase dt and keep the evolution and limit scaling factor to 5.0
  115. dt *= static_cast<value_type>(9)/static_cast<value_type>(10) *
  116. pow(error, static_cast<value_type>(-1) / stepper_order);
  117. if(m_max_dt != static_cast<time_type >(0))
  118. // limit to maximal stepsize
  119. dt = detail::min_abs(dt, m_max_dt);
  120. }
  121. return dt;
  122. }
  123. bool check_step_size_limit(const time_type dt)
  124. {
  125. if(m_max_dt != static_cast<time_type >(0))
  126. return detail::less_eq_with_sign(dt, m_max_dt, dt);
  127. return true;
  128. }
  129. time_type get_max_dt() { return m_max_dt; }
  130. protected:
  131. time_type m_max_dt;
  132. };
  133. /*
  134. * error stepper category dispatcher
  135. */
  136. template<
  137. class ErrorStepper ,
  138. class ErrorChecker = default_error_checker< typename ErrorStepper::value_type ,
  139. typename ErrorStepper::algebra_type ,
  140. typename ErrorStepper::operations_type > ,
  141. class StepAdjuster = default_step_adjuster< typename ErrorStepper::value_type ,
  142. typename ErrorStepper::time_type > ,
  143. class Resizer = typename ErrorStepper::resizer_type ,
  144. class ErrorStepperCategory = typename ErrorStepper::stepper_category
  145. >
  146. class controlled_runge_kutta ;
  147. /*
  148. * explicit stepper version
  149. *
  150. * this class introduces the following try_step overloads
  151. * try_step( sys , x , t , dt )
  152. * try_step( sys , x , dxdt , t , dt )
  153. * try_step( sys , in , t , out , dt )
  154. * try_step( sys , in , dxdt , t , out , dt )
  155. */
  156. /**
  157. * \brief Implements step size control for Runge-Kutta steppers with error
  158. * estimation.
  159. *
  160. * This class implements the step size control for standard Runge-Kutta
  161. * steppers with error estimation.
  162. *
  163. * \tparam ErrorStepper The stepper type with error estimation, has to fulfill the ErrorStepper concept.
  164. * \tparam ErrorChecker The error checker
  165. * \tparam Resizer The resizer policy type.
  166. */
  167. template<
  168. class ErrorStepper,
  169. class ErrorChecker,
  170. class StepAdjuster,
  171. class Resizer
  172. >
  173. class controlled_runge_kutta< ErrorStepper , ErrorChecker , StepAdjuster, Resizer ,
  174. explicit_error_stepper_tag >
  175. {
  176. public:
  177. typedef ErrorStepper stepper_type;
  178. typedef typename stepper_type::state_type state_type;
  179. typedef typename stepper_type::value_type value_type;
  180. typedef typename stepper_type::deriv_type deriv_type;
  181. typedef typename stepper_type::time_type time_type;
  182. typedef typename stepper_type::algebra_type algebra_type;
  183. typedef typename stepper_type::operations_type operations_type;
  184. typedef Resizer resizer_type;
  185. typedef ErrorChecker error_checker_type;
  186. typedef StepAdjuster step_adjuster_type;
  187. typedef explicit_controlled_stepper_tag stepper_category;
  188. #ifndef DOXYGEN_SKIP
  189. typedef typename stepper_type::wrapped_state_type wrapped_state_type;
  190. typedef typename stepper_type::wrapped_deriv_type wrapped_deriv_type;
  191. typedef controlled_runge_kutta< ErrorStepper , ErrorChecker , StepAdjuster ,
  192. Resizer , explicit_error_stepper_tag > controlled_stepper_type;
  193. #endif //DOXYGEN_SKIP
  194. /**
  195. * \brief Constructs the controlled Runge-Kutta stepper.
  196. * \param error_checker An instance of the error checker.
  197. * \param stepper An instance of the underlying stepper.
  198. */
  199. controlled_runge_kutta(
  200. const error_checker_type &error_checker = error_checker_type( ) ,
  201. const step_adjuster_type &step_adjuster = step_adjuster_type() ,
  202. const stepper_type &stepper = stepper_type( )
  203. )
  204. : m_stepper(stepper), m_error_checker(error_checker) , m_step_adjuster(step_adjuster)
  205. { }
  206. /*
  207. * Version 1 : try_step( sys , x , t , dt )
  208. *
  209. * The overloads are needed to solve the forwarding problem
  210. */
  211. /**
  212. * \brief Tries to perform one step.
  213. *
  214. * This method tries to do one step with step size dt. If the error estimate
  215. * is to large, the step is rejected and the method returns fail and the
  216. * step size dt is reduced. If the error estimate is acceptably small, the
  217. * step is performed, success is returned and dt might be increased to make
  218. * the steps as large as possible. This method also updates t if a step is
  219. * performed.
  220. *
  221. * \param system The system function to solve, hence the r.h.s. of the ODE. It must fulfill the
  222. * Simple System concept.
  223. * \param x The state of the ODE which should be solved. Overwritten if
  224. * the step is successful.
  225. * \param t The value of the time. Updated if the step is successful.
  226. * \param dt The step size. Updated.
  227. * \return success if the step was accepted, fail otherwise.
  228. */
  229. template< class System , class StateInOut >
  230. controlled_step_result try_step( System system , StateInOut &x , time_type &t , time_type &dt )
  231. {
  232. return try_step_v1( system , x , t, dt );
  233. }
  234. /**
  235. * \brief Tries to perform one step. Solves the forwarding problem and
  236. * allows for using boost range as state_type.
  237. *
  238. * This method tries to do one step with step size dt. If the error estimate
  239. * is to large, the step is rejected and the method returns fail and the
  240. * step size dt is reduced. If the error estimate is acceptably small, the
  241. * step is performed, success is returned and dt might be increased to make
  242. * the steps as large as possible. This method also updates t if a step is
  243. * performed.
  244. *
  245. * \param system The system function to solve, hence the r.h.s. of the ODE. It must fulfill the
  246. * Simple System concept.
  247. * \param x The state of the ODE which should be solved. Overwritten if
  248. * the step is successful. Can be a boost range.
  249. * \param t The value of the time. Updated if the step is successful.
  250. * \param dt The step size. Updated.
  251. * \return success if the step was accepted, fail otherwise.
  252. */
  253. template< class System , class StateInOut >
  254. controlled_step_result try_step( System system , const StateInOut &x , time_type &t , time_type &dt )
  255. {
  256. return try_step_v1( system , x , t, dt );
  257. }
  258. /*
  259. * Version 2 : try_step( sys , x , dxdt , t , dt )
  260. *
  261. * this version does not solve the forwarding problem, boost.range can not be used
  262. */
  263. /**
  264. * \brief Tries to perform one step.
  265. *
  266. * This method tries to do one step with step size dt. If the error estimate
  267. * is to large, the step is rejected and the method returns fail and the
  268. * step size dt is reduced. If the error estimate is acceptably small, the
  269. * step is performed, success is returned and dt might be increased to make
  270. * the steps as large as possible. This method also updates t if a step is
  271. * performed.
  272. *
  273. * \param system The system function to solve, hence the r.h.s. of the ODE. It must fulfill the
  274. * Simple System concept.
  275. * \param x The state of the ODE which should be solved. Overwritten if
  276. * the step is successful.
  277. * \param dxdt The derivative of state.
  278. * \param t The value of the time. Updated if the step is successful.
  279. * \param dt The step size. Updated.
  280. * \return success if the step was accepted, fail otherwise.
  281. */
  282. template< class System , class StateInOut , class DerivIn >
  283. controlled_step_result try_step( System system , StateInOut &x , const DerivIn &dxdt , time_type &t , time_type &dt )
  284. {
  285. m_xnew_resizer.adjust_size( x , detail::bind( &controlled_runge_kutta::template resize_m_xnew_impl< StateInOut > , detail::ref( *this ) , detail::_1 ) );
  286. controlled_step_result res = try_step( system , x , dxdt , t , m_xnew.m_v , dt );
  287. if( res == success )
  288. {
  289. boost::numeric::odeint::copy( m_xnew.m_v , x );
  290. }
  291. return res;
  292. }
  293. /*
  294. * Version 3 : try_step( sys , in , t , out , dt )
  295. *
  296. * this version does not solve the forwarding problem, boost.range can not be used
  297. *
  298. * the disable is needed to avoid ambiguous overloads if state_type = time_type
  299. */
  300. /**
  301. * \brief Tries to perform one step.
  302. *
  303. * \note This method is disabled if state_type=time_type to avoid ambiguity.
  304. *
  305. * This method tries to do one step with step size dt. If the error estimate
  306. * is to large, the step is rejected and the method returns fail and the
  307. * step size dt is reduced. If the error estimate is acceptably small, the
  308. * step is performed, success is returned and dt might be increased to make
  309. * the steps as large as possible. This method also updates t if a step is
  310. * performed.
  311. *
  312. * \param system The system function to solve, hence the r.h.s. of the ODE. It must fulfill the
  313. * Simple System concept.
  314. * \param in The state of the ODE which should be solved.
  315. * \param t The value of the time. Updated if the step is successful.
  316. * \param out Used to store the result of the step.
  317. * \param dt The step size. Updated.
  318. * \return success if the step was accepted, fail otherwise.
  319. */
  320. template< class System , class StateIn , class StateOut >
  321. typename boost::disable_if< boost::is_same< StateIn , time_type > , controlled_step_result >::type
  322. try_step( System system , const StateIn &in , time_type &t , StateOut &out , time_type &dt )
  323. {
  324. typename odeint::unwrap_reference< System >::type &sys = system;
  325. m_dxdt_resizer.adjust_size( in , detail::bind( &controlled_runge_kutta::template resize_m_dxdt_impl< StateIn > , detail::ref( *this ) , detail::_1 ) );
  326. sys( in , m_dxdt.m_v , t );
  327. return try_step( system , in , m_dxdt.m_v , t , out , dt );
  328. }
  329. /*
  330. * Version 4 : try_step( sys , in , dxdt , t , out , dt )
  331. *
  332. * this version does not solve the forwarding problem, boost.range can not be used
  333. */
  334. /**
  335. * \brief Tries to perform one step.
  336. *
  337. * This method tries to do one step with step size dt. If the error estimate
  338. * is to large, the step is rejected and the method returns fail and the
  339. * step size dt is reduced. If the error estimate is acceptably small, the
  340. * step is performed, success is returned and dt might be increased to make
  341. * the steps as large as possible. This method also updates t if a step is
  342. * performed.
  343. *
  344. * \param system The system function to solve, hence the r.h.s. of the ODE. It must fulfill the
  345. * Simple System concept.
  346. * \param in The state of the ODE which should be solved.
  347. * \param dxdt The derivative of state.
  348. * \param t The value of the time. Updated if the step is successful.
  349. * \param out Used to store the result of the step.
  350. * \param dt The step size. Updated.
  351. * \return success if the step was accepted, fail otherwise.
  352. */
  353. template< class System , class StateIn , class DerivIn , class StateOut >
  354. controlled_step_result try_step( System system , const StateIn &in , const DerivIn &dxdt , time_type &t , StateOut &out , time_type &dt )
  355. {
  356. unwrapped_step_adjuster &step_adjuster = m_step_adjuster;
  357. if( !step_adjuster.check_step_size_limit(dt) )
  358. {
  359. // given dt was above step size limit - adjust and return fail;
  360. dt = step_adjuster.get_max_dt();
  361. return fail;
  362. }
  363. m_xerr_resizer.adjust_size( in , detail::bind( &controlled_runge_kutta::template resize_m_xerr_impl< StateIn > , detail::ref( *this ) , detail::_1 ) );
  364. // do one step with error calculation
  365. m_stepper.do_step( system , in , dxdt , t , out , dt , m_xerr.m_v );
  366. value_type max_rel_err = m_error_checker.error( m_stepper.algebra() , in , dxdt , m_xerr.m_v , dt );
  367. if( max_rel_err > 1.0 )
  368. {
  369. // error too big, decrease step size and reject this step
  370. dt = step_adjuster.decrease_step(dt, max_rel_err, m_stepper.error_order());
  371. return fail;
  372. } else
  373. {
  374. // otherwise, increase step size and accept
  375. t += dt;
  376. dt = step_adjuster.increase_step(dt, max_rel_err, m_stepper.stepper_order());
  377. return success;
  378. }
  379. }
  380. /**
  381. * \brief Adjust the size of all temporaries in the stepper manually.
  382. * \param x A state from which the size of the temporaries to be resized is deduced.
  383. */
  384. template< class StateType >
  385. void adjust_size( const StateType &x )
  386. {
  387. resize_m_xerr_impl( x );
  388. resize_m_dxdt_impl( x );
  389. resize_m_xnew_impl( x );
  390. m_stepper.adjust_size( x );
  391. }
  392. /**
  393. * \brief Returns the instance of the underlying stepper.
  394. * \returns The instance of the underlying stepper.
  395. */
  396. stepper_type& stepper( void )
  397. {
  398. return m_stepper;
  399. }
  400. /**
  401. * \brief Returns the instance of the underlying stepper.
  402. * \returns The instance of the underlying stepper.
  403. */
  404. const stepper_type& stepper( void ) const
  405. {
  406. return m_stepper;
  407. }
  408. private:
  409. template< class System , class StateInOut >
  410. controlled_step_result try_step_v1( System system , StateInOut &x , time_type &t , time_type &dt )
  411. {
  412. typename odeint::unwrap_reference< System >::type &sys = system;
  413. m_dxdt_resizer.adjust_size( x , detail::bind( &controlled_runge_kutta::template resize_m_dxdt_impl< StateInOut > , detail::ref( *this ) , detail::_1 ) );
  414. sys( x , m_dxdt.m_v ,t );
  415. return try_step( system , x , m_dxdt.m_v , t , dt );
  416. }
  417. template< class StateIn >
  418. bool resize_m_xerr_impl( const StateIn &x )
  419. {
  420. return adjust_size_by_resizeability( m_xerr , x , typename is_resizeable<state_type>::type() );
  421. }
  422. template< class StateIn >
  423. bool resize_m_dxdt_impl( const StateIn &x )
  424. {
  425. return adjust_size_by_resizeability( m_dxdt , x , typename is_resizeable<deriv_type>::type() );
  426. }
  427. template< class StateIn >
  428. bool resize_m_xnew_impl( const StateIn &x )
  429. {
  430. return adjust_size_by_resizeability( m_xnew , x , typename is_resizeable<state_type>::type() );
  431. }
  432. stepper_type m_stepper;
  433. error_checker_type m_error_checker;
  434. step_adjuster_type m_step_adjuster;
  435. typedef typename unwrap_reference< step_adjuster_type >::type unwrapped_step_adjuster;
  436. resizer_type m_dxdt_resizer;
  437. resizer_type m_xerr_resizer;
  438. resizer_type m_xnew_resizer;
  439. wrapped_deriv_type m_dxdt;
  440. wrapped_state_type m_xerr;
  441. wrapped_state_type m_xnew;
  442. };
  443. /*
  444. * explicit stepper fsal version
  445. *
  446. * the class introduces the following try_step overloads
  447. * try_step( sys , x , t , dt )
  448. * try_step( sys , in , t , out , dt )
  449. * try_step( sys , x , dxdt , t , dt )
  450. * try_step( sys , in , dxdt_in , t , out , dxdt_out , dt )
  451. */
  452. /**
  453. * \brief Implements step size control for Runge-Kutta FSAL steppers with
  454. * error estimation.
  455. *
  456. * This class implements the step size control for FSAL Runge-Kutta
  457. * steppers with error estimation.
  458. *
  459. * \tparam ErrorStepper The stepper type with error estimation, has to fulfill the ErrorStepper concept.
  460. * \tparam ErrorChecker The error checker
  461. * \tparam Resizer The resizer policy type.
  462. */
  463. template<
  464. class ErrorStepper ,
  465. class ErrorChecker ,
  466. class StepAdjuster ,
  467. class Resizer
  468. >
  469. class controlled_runge_kutta< ErrorStepper , ErrorChecker , StepAdjuster , Resizer , explicit_error_stepper_fsal_tag >
  470. {
  471. public:
  472. typedef ErrorStepper stepper_type;
  473. typedef typename stepper_type::state_type state_type;
  474. typedef typename stepper_type::value_type value_type;
  475. typedef typename stepper_type::deriv_type deriv_type;
  476. typedef typename stepper_type::time_type time_type;
  477. typedef typename stepper_type::algebra_type algebra_type;
  478. typedef typename stepper_type::operations_type operations_type;
  479. typedef Resizer resizer_type;
  480. typedef ErrorChecker error_checker_type;
  481. typedef StepAdjuster step_adjuster_type;
  482. typedef explicit_controlled_stepper_fsal_tag stepper_category;
  483. #ifndef DOXYGEN_SKIP
  484. typedef typename stepper_type::wrapped_state_type wrapped_state_type;
  485. typedef typename stepper_type::wrapped_deriv_type wrapped_deriv_type;
  486. typedef controlled_runge_kutta< ErrorStepper , ErrorChecker , StepAdjuster , Resizer , explicit_error_stepper_tag > controlled_stepper_type;
  487. #endif // DOXYGEN_SKIP
  488. /**
  489. * \brief Constructs the controlled Runge-Kutta stepper.
  490. * \param error_checker An instance of the error checker.
  491. * \param stepper An instance of the underlying stepper.
  492. */
  493. controlled_runge_kutta(
  494. const error_checker_type &error_checker = error_checker_type() ,
  495. const step_adjuster_type &step_adjuster = step_adjuster_type() ,
  496. const stepper_type &stepper = stepper_type()
  497. )
  498. : m_stepper( stepper ) , m_error_checker( error_checker ) , m_step_adjuster(step_adjuster) ,
  499. m_first_call( true )
  500. { }
  501. /*
  502. * Version 1 : try_step( sys , x , t , dt )
  503. *
  504. * The two overloads are needed in order to solve the forwarding problem
  505. */
  506. /**
  507. * \brief Tries to perform one step.
  508. *
  509. * This method tries to do one step with step size dt. If the error estimate
  510. * is to large, the step is rejected and the method returns fail and the
  511. * step size dt is reduced. If the error estimate is acceptably small, the
  512. * step is performed, success is returned and dt might be increased to make
  513. * the steps as large as possible. This method also updates t if a step is
  514. * performed.
  515. *
  516. * \param system The system function to solve, hence the r.h.s. of the ODE. It must fulfill the
  517. * Simple System concept.
  518. * \param x The state of the ODE which should be solved. Overwritten if
  519. * the step is successful.
  520. * \param t The value of the time. Updated if the step is successful.
  521. * \param dt The step size. Updated.
  522. * \return success if the step was accepted, fail otherwise.
  523. */
  524. template< class System , class StateInOut >
  525. controlled_step_result try_step( System system , StateInOut &x , time_type &t , time_type &dt )
  526. {
  527. return try_step_v1( system , x , t , dt );
  528. }
  529. /**
  530. * \brief Tries to perform one step. Solves the forwarding problem and
  531. * allows for using boost range as state_type.
  532. *
  533. * This method tries to do one step with step size dt. If the error estimate
  534. * is to large, the step is rejected and the method returns fail and the
  535. * step size dt is reduced. If the error estimate is acceptably small, the
  536. * step is performed, success is returned and dt might be increased to make
  537. * the steps as large as possible. This method also updates t if a step is
  538. * performed.
  539. *
  540. * \param system The system function to solve, hence the r.h.s. of the ODE. It must fulfill the
  541. * Simple System concept.
  542. * \param x The state of the ODE which should be solved. Overwritten if
  543. * the step is successful. Can be a boost range.
  544. * \param t The value of the time. Updated if the step is successful.
  545. * \param dt The step size. Updated.
  546. * \return success if the step was accepted, fail otherwise.
  547. */
  548. template< class System , class StateInOut >
  549. controlled_step_result try_step( System system , const StateInOut &x , time_type &t , time_type &dt )
  550. {
  551. return try_step_v1( system , x , t , dt );
  552. }
  553. /*
  554. * Version 2 : try_step( sys , in , t , out , dt );
  555. *
  556. * This version does not solve the forwarding problem, boost::range can not be used.
  557. *
  558. * The disabler is needed to solve ambiguous overloads
  559. */
  560. /**
  561. * \brief Tries to perform one step.
  562. *
  563. * \note This method is disabled if state_type=time_type to avoid ambiguity.
  564. *
  565. * This method tries to do one step with step size dt. If the error estimate
  566. * is to large, the step is rejected and the method returns fail and the
  567. * step size dt is reduced. If the error estimate is acceptably small, the
  568. * step is performed, success is returned and dt might be increased to make
  569. * the steps as large as possible. This method also updates t if a step is
  570. * performed.
  571. *
  572. * \param system The system function to solve, hence the r.h.s. of the ODE. It must fulfill the
  573. * Simple System concept.
  574. * \param in The state of the ODE which should be solved.
  575. * \param t The value of the time. Updated if the step is successful.
  576. * \param out Used to store the result of the step.
  577. * \param dt The step size. Updated.
  578. * \return success if the step was accepted, fail otherwise.
  579. */
  580. template< class System , class StateIn , class StateOut >
  581. typename boost::disable_if< boost::is_same< StateIn , time_type > , controlled_step_result >::type
  582. try_step( System system , const StateIn &in , time_type &t , StateOut &out , time_type &dt )
  583. {
  584. if( m_dxdt_resizer.adjust_size( in , detail::bind( &controlled_runge_kutta::template resize_m_dxdt_impl< StateIn > , detail::ref( *this ) , detail::_1 ) ) || m_first_call )
  585. {
  586. initialize( system , in , t );
  587. }
  588. return try_step( system , in , m_dxdt.m_v , t , out , dt );
  589. }
  590. /*
  591. * Version 3 : try_step( sys , x , dxdt , t , dt )
  592. *
  593. * This version does not solve the forwarding problem, boost::range can not be used.
  594. */
  595. /**
  596. * \brief Tries to perform one step.
  597. *
  598. * This method tries to do one step with step size dt. If the error estimate
  599. * is to large, the step is rejected and the method returns fail and the
  600. * step size dt is reduced. If the error estimate is acceptably small, the
  601. * step is performed, success is returned and dt might be increased to make
  602. * the steps as large as possible. This method also updates t if a step is
  603. * performed.
  604. *
  605. * \param system The system function to solve, hence the r.h.s. of the ODE. It must fulfill the
  606. * Simple System concept.
  607. * \param x The state of the ODE which should be solved. Overwritten if
  608. * the step is successful.
  609. * \param dxdt The derivative of state.
  610. * \param t The value of the time. Updated if the step is successful.
  611. * \param dt The step size. Updated.
  612. * \return success if the step was accepted, fail otherwise.
  613. */
  614. template< class System , class StateInOut , class DerivInOut >
  615. controlled_step_result try_step( System system , StateInOut &x , DerivInOut &dxdt , time_type &t , time_type &dt )
  616. {
  617. m_xnew_resizer.adjust_size( x , detail::bind( &controlled_runge_kutta::template resize_m_xnew_impl< StateInOut > , detail::ref( *this ) , detail::_1 ) );
  618. m_dxdt_new_resizer.adjust_size( x , detail::bind( &controlled_runge_kutta::template resize_m_dxdt_new_impl< StateInOut > , detail::ref( *this ) , detail::_1 ) );
  619. controlled_step_result res = try_step( system , x , dxdt , t , m_xnew.m_v , m_dxdtnew.m_v , dt );
  620. if( res == success )
  621. {
  622. boost::numeric::odeint::copy( m_xnew.m_v , x );
  623. boost::numeric::odeint::copy( m_dxdtnew.m_v , dxdt );
  624. }
  625. return res;
  626. }
  627. /*
  628. * Version 4 : try_step( sys , in , dxdt_in , t , out , dxdt_out , dt )
  629. *
  630. * This version does not solve the forwarding problem, boost::range can not be used.
  631. */
  632. /**
  633. * \brief Tries to perform one step.
  634. *
  635. * This method tries to do one step with step size dt. If the error estimate
  636. * is to large, the step is rejected and the method returns fail and the
  637. * step size dt is reduced. If the error estimate is acceptably small, the
  638. * step is performed, success is returned and dt might be increased to make
  639. * the steps as large as possible. This method also updates t if a step is
  640. * performed.
  641. *
  642. * \param system The system function to solve, hence the r.h.s. of the ODE. It must fulfill the
  643. * Simple System concept.
  644. * \param in The state of the ODE which should be solved.
  645. * \param dxdt The derivative of state.
  646. * \param t The value of the time. Updated if the step is successful.
  647. * \param out Used to store the result of the step.
  648. * \param dt The step size. Updated.
  649. * \return success if the step was accepted, fail otherwise.
  650. */
  651. template< class System , class StateIn , class DerivIn , class StateOut , class DerivOut >
  652. controlled_step_result try_step( System system , const StateIn &in , const DerivIn &dxdt_in , time_type &t ,
  653. StateOut &out , DerivOut &dxdt_out , time_type &dt )
  654. {
  655. unwrapped_step_adjuster &step_adjuster = m_step_adjuster;
  656. if( !step_adjuster.check_step_size_limit(dt) )
  657. {
  658. // given dt was above step size limit - adjust and return fail;
  659. dt = step_adjuster.get_max_dt();
  660. return fail;
  661. }
  662. m_xerr_resizer.adjust_size( in , detail::bind( &controlled_runge_kutta::template resize_m_xerr_impl< StateIn > , detail::ref( *this ) , detail::_1 ) );
  663. //fsal: m_stepper.get_dxdt( dxdt );
  664. //fsal: m_stepper.do_step( sys , x , dxdt , t , dt , m_x_err );
  665. m_stepper.do_step( system , in , dxdt_in , t , out , dxdt_out , dt , m_xerr.m_v );
  666. // this potentially overwrites m_x_err! (standard_error_checker does, at least)
  667. value_type max_rel_err = m_error_checker.error( m_stepper.algebra() , in , dxdt_in , m_xerr.m_v , dt );
  668. if( max_rel_err > 1.0 )
  669. {
  670. // error too big, decrease step size and reject this step
  671. dt = step_adjuster.decrease_step(dt, max_rel_err, m_stepper.error_order());
  672. return fail;
  673. }
  674. // otherwise, increase step size and accept
  675. t += dt;
  676. dt = step_adjuster.increase_step(dt, max_rel_err, m_stepper.stepper_order());
  677. return success;
  678. }
  679. /**
  680. * \brief Resets the internal state of the underlying FSAL stepper.
  681. */
  682. void reset( void )
  683. {
  684. m_first_call = true;
  685. }
  686. /**
  687. * \brief Initializes the internal state storing an internal copy of the derivative.
  688. *
  689. * \param deriv The initial derivative of the ODE.
  690. */
  691. template< class DerivIn >
  692. void initialize( const DerivIn &deriv )
  693. {
  694. boost::numeric::odeint::copy( deriv , m_dxdt.m_v );
  695. m_first_call = false;
  696. }
  697. /**
  698. * \brief Initializes the internal state storing an internal copy of the derivative.
  699. *
  700. * \param system The system function to solve, hence the r.h.s. of the ODE. It must fulfill the
  701. * Simple System concept.
  702. * \param x The initial state of the ODE which should be solved.
  703. * \param t The initial time.
  704. */
  705. template< class System , class StateIn >
  706. void initialize( System system , const StateIn &x , time_type t )
  707. {
  708. typename odeint::unwrap_reference< System >::type &sys = system;
  709. sys( x , m_dxdt.m_v , t );
  710. m_first_call = false;
  711. }
  712. /**
  713. * \brief Returns true if the stepper has been initialized, false otherwise.
  714. *
  715. * \return true, if the stepper has been initialized, false otherwise.
  716. */
  717. bool is_initialized( void ) const
  718. {
  719. return ! m_first_call;
  720. }
  721. /**
  722. * \brief Adjust the size of all temporaries in the stepper manually.
  723. * \param x A state from which the size of the temporaries to be resized is deduced.
  724. */
  725. template< class StateType >
  726. void adjust_size( const StateType &x )
  727. {
  728. resize_m_xerr_impl( x );
  729. resize_m_dxdt_impl( x );
  730. resize_m_dxdt_new_impl( x );
  731. resize_m_xnew_impl( x );
  732. }
  733. /**
  734. * \brief Returns the instance of the underlying stepper.
  735. * \returns The instance of the underlying stepper.
  736. */
  737. stepper_type& stepper( void )
  738. {
  739. return m_stepper;
  740. }
  741. /**
  742. * \brief Returns the instance of the underlying stepper.
  743. * \returns The instance of the underlying stepper.
  744. */
  745. const stepper_type& stepper( void ) const
  746. {
  747. return m_stepper;
  748. }
  749. private:
  750. template< class StateIn >
  751. bool resize_m_xerr_impl( const StateIn &x )
  752. {
  753. return adjust_size_by_resizeability( m_xerr , x , typename is_resizeable<state_type>::type() );
  754. }
  755. template< class StateIn >
  756. bool resize_m_dxdt_impl( const StateIn &x )
  757. {
  758. return adjust_size_by_resizeability( m_dxdt , x , typename is_resizeable<deriv_type>::type() );
  759. }
  760. template< class StateIn >
  761. bool resize_m_dxdt_new_impl( const StateIn &x )
  762. {
  763. return adjust_size_by_resizeability( m_dxdtnew , x , typename is_resizeable<deriv_type>::type() );
  764. }
  765. template< class StateIn >
  766. bool resize_m_xnew_impl( const StateIn &x )
  767. {
  768. return adjust_size_by_resizeability( m_xnew , x , typename is_resizeable<state_type>::type() );
  769. }
  770. template< class System , class StateInOut >
  771. controlled_step_result try_step_v1( System system , StateInOut &x , time_type &t , time_type &dt )
  772. {
  773. if( m_dxdt_resizer.adjust_size( x , detail::bind( &controlled_runge_kutta::template resize_m_dxdt_impl< StateInOut > , detail::ref( *this ) , detail::_1 ) ) || m_first_call )
  774. {
  775. initialize( system , x , t );
  776. }
  777. return try_step( system , x , m_dxdt.m_v , t , dt );
  778. }
  779. stepper_type m_stepper;
  780. error_checker_type m_error_checker;
  781. step_adjuster_type m_step_adjuster;
  782. typedef typename unwrap_reference< step_adjuster_type >::type unwrapped_step_adjuster;
  783. resizer_type m_dxdt_resizer;
  784. resizer_type m_xerr_resizer;
  785. resizer_type m_xnew_resizer;
  786. resizer_type m_dxdt_new_resizer;
  787. wrapped_deriv_type m_dxdt;
  788. wrapped_state_type m_xerr;
  789. wrapped_state_type m_xnew;
  790. wrapped_deriv_type m_dxdtnew;
  791. bool m_first_call;
  792. };
  793. /********** DOXYGEN **********/
  794. /**** DEFAULT ERROR CHECKER ****/
  795. /**
  796. * \class default_error_checker
  797. * \brief The default error checker to be used with Runge-Kutta error steppers
  798. *
  799. * This class provides the default mechanism to compare the error estimates
  800. * reported by Runge-Kutta error steppers with user defined error bounds.
  801. * It is used by the controlled_runge_kutta steppers.
  802. *
  803. * \tparam Value The value type.
  804. * \tparam Time The time type.
  805. * \tparam Algebra The algebra type.
  806. * \tparam Operations The operations type.
  807. */
  808. /**
  809. * \fn default_error_checker( value_type eps_abs , value_type eps_rel , value_type a_x , value_type a_dxdt ,
  810. * time_type max_dt)
  811. * \brief Constructs the error checker.
  812. *
  813. * The error is calculated as follows: ????
  814. *
  815. * \param eps_abs Absolute tolerance level.
  816. * \param eps_rel Relative tolerance level.
  817. * \param a_x Factor for the weight of the state.
  818. * \param a_dxdt Factor for the weight of the derivative.
  819. * \param max_dt Maximum allowed step size.
  820. */
  821. /**
  822. * \fn error( const State &x_old , const Deriv &dxdt_old , Err &x_err , time_type dt ) const
  823. * \brief Calculates the error level.
  824. *
  825. * If the returned error level is greater than 1, the estimated error was
  826. * larger than the permitted error bounds and the step should be repeated
  827. * with a smaller step size.
  828. *
  829. * \param x_old State at the beginning of the step.
  830. * \param dxdt_old Derivative at the beginning of the step.
  831. * \param x_err Error estimate.
  832. * \param dt Time step.
  833. * \return error
  834. */
  835. /**
  836. * \fn error( algebra_type &algebra , const State &x_old , const Deriv &dxdt_old , Err &x_err , time_type dt ) const
  837. * \brief Calculates the error level using a given algebra.
  838. *
  839. * If the returned error level is greater than 1, the estimated error was
  840. * larger than the permitted error bounds and the step should be repeated
  841. * with a smaller step size.
  842. *
  843. * \param algebra The algebra used for calculation of the error.
  844. * \param x_old State at the beginning of the step.
  845. * \param dxdt_old Derivative at the beginning of the step.
  846. * \param x_err Error estimate.
  847. * \param dt Time step.
  848. * \return error
  849. */
  850. /**
  851. * \fn time_type decrease_step(const time_type dt, const value_type error, const int error_order)
  852. * \brief Returns a decreased step size based on the given error and order
  853. *
  854. * Calculates a new smaller step size based on the given error and its order.
  855. *
  856. * \param dt The old step size.
  857. * \param error The computed error estimate.
  858. * \param error_order The error order of the stepper.
  859. * \return dt_new The new, reduced step size.
  860. */
  861. /**
  862. * \fn time_type increase_step(const time_type dt, const value_type error, const int error_order)
  863. * \brief Returns an increased step size based on the given error and order.
  864. *
  865. * Calculates a new bigger step size based on the given error and its order. If max_dt != 0, the
  866. * new step size is limited to max_dt.
  867. *
  868. * \param dt The old step size.
  869. * \param error The computed error estimate.
  870. * \param error_order The order of the stepper.
  871. * \return dt_new The new, increased step size.
  872. */
  873. } // odeint
  874. } // numeric
  875. } // boost
  876. #endif // BOOST_NUMERIC_ODEINT_STEPPER_CONTROLLED_RUNGE_KUTTA_HPP_INCLUDED