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