adams_moulton.hpp 5.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195
  1. /*
  2. [auto_generated]
  3. boost/numeric/odeint/stepper/adams_moulton.hpp
  4. [begin_description]
  5. Implementation of the Adams-Moulton method. This is method is not a real stepper, it is more a helper class
  6. which computes the corrector step in the Adams-Bashforth-Moulton method.
  7. [end_description]
  8. Copyright 2009-2011 Karsten Ahnert
  9. Copyright 2009-2011 Mario Mulansky
  10. Distributed under the Boost Software License, Version 1.0.
  11. (See accompanying file LICENSE_1_0.txt or
  12. copy at http://www.boost.org/LICENSE_1_0.txt)
  13. */
  14. #ifndef BOOST_NUMERIC_ODEINT_STEPPER_ADAMS_MOULTON_HPP_INCLUDED
  15. #define BOOST_NUMERIC_ODEINT_STEPPER_ADAMS_MOULTON_HPP_INCLUDED
  16. #include <boost/numeric/odeint/util/bind.hpp>
  17. #include <boost/numeric/odeint/algebra/range_algebra.hpp>
  18. #include <boost/numeric/odeint/algebra/default_operations.hpp>
  19. #include <boost/numeric/odeint/util/state_wrapper.hpp>
  20. #include <boost/numeric/odeint/util/is_resizeable.hpp>
  21. #include <boost/numeric/odeint/util/resizer.hpp>
  22. #include <boost/numeric/odeint/stepper/stepper_categories.hpp>
  23. #include <boost/numeric/odeint/stepper/runge_kutta4_classic.hpp>
  24. #include <boost/numeric/odeint/stepper/detail/adams_moulton_call_algebra.hpp>
  25. #include <boost/numeric/odeint/stepper/detail/adams_moulton_coefficients.hpp>
  26. #include <boost/numeric/odeint/stepper/detail/rotating_buffer.hpp>
  27. namespace boost {
  28. namespace numeric {
  29. namespace odeint {
  30. /*
  31. * Static implicit Adams-Moulton multistep-solver without step size control and without dense output.
  32. */
  33. template<
  34. size_t Steps ,
  35. class State ,
  36. class Value = double ,
  37. class Deriv = State ,
  38. class Time = Value ,
  39. class Algebra = range_algebra ,
  40. class Operations = default_operations ,
  41. class Resizer = initially_resizer
  42. >
  43. class adams_moulton
  44. {
  45. private:
  46. public :
  47. typedef State state_type;
  48. typedef state_wrapper< state_type > wrapped_state_type;
  49. typedef Value value_type;
  50. typedef Deriv deriv_type;
  51. typedef state_wrapper< deriv_type > wrapped_deriv_type;
  52. typedef Time time_type;
  53. typedef Algebra algebra_type;
  54. typedef Operations operations_type;
  55. typedef Resizer resizer_type;
  56. typedef stepper_tag stepper_category;
  57. typedef adams_moulton< Steps , State , Value , Deriv , Time , Algebra , Operations , Resizer > stepper_type;
  58. static const size_t steps = Steps;
  59. typedef unsigned short order_type;
  60. static const order_type order_value = steps + 1;
  61. typedef detail::rotating_buffer< wrapped_deriv_type , steps > step_storage_type;
  62. adams_moulton( )
  63. : m_coefficients() , m_dxdt() , m_resizer() ,
  64. m_algebra_instance() , m_algebra( m_algebra_instance )
  65. { }
  66. adams_moulton( algebra_type &algebra )
  67. : m_coefficients() , m_dxdt() , m_resizer() ,
  68. m_algebra_instance() , m_algebra( algebra )
  69. { }
  70. adams_moulton& operator=( const adams_moulton &stepper )
  71. {
  72. m_dxdt = stepper.m_dxdt;
  73. m_resizer = stepper.m_resizer;
  74. m_algebra = stepper.m_algebra;
  75. return *this;
  76. }
  77. order_type order( void ) const { return order_value; }
  78. /*
  79. * Version 1 : do_step( system , x , t , dt , buf );
  80. *
  81. * solves the forwarding problem
  82. */
  83. template< class System , class StateInOut , class ABBuf >
  84. void do_step( System system , StateInOut &in , time_type t , time_type dt , const ABBuf &buf )
  85. {
  86. do_step( system , in , t , in , dt , buf );
  87. }
  88. template< class System , class StateInOut , class ABBuf >
  89. void do_step( System system , const StateInOut &in , time_type t , time_type dt , const ABBuf &buf )
  90. {
  91. do_step( system , in , t , in , dt , buf );
  92. }
  93. /*
  94. * Version 2 : do_step( system , in , t , out , dt , buf );
  95. *
  96. * solves the forwarding problem
  97. */
  98. template< class System , class StateIn , class StateOut , class ABBuf >
  99. void do_step( System system , const StateIn &in , time_type t , StateOut &out , time_type dt , const ABBuf &buf )
  100. {
  101. typename odeint::unwrap_reference< System >::type &sys = system;
  102. m_resizer.adjust_size( in , detail::bind( &stepper_type::template resize_impl<StateIn> , detail::ref( *this ) , detail::_1 ) );
  103. sys( in , m_dxdt.m_v , t );
  104. detail::adams_moulton_call_algebra< steps , algebra_type , operations_type >()( m_algebra , in , out , m_dxdt.m_v , buf , m_coefficients , dt );
  105. }
  106. template< class System , class StateIn , class StateOut , class ABBuf >
  107. void do_step( System system , const StateIn &in , time_type t , const StateOut &out , time_type dt , const ABBuf &buf )
  108. {
  109. typename odeint::unwrap_reference< System >::type &sys = system;
  110. m_resizer.adjust_size( in , detail::bind( &stepper_type::template resize_impl<StateIn> , detail::ref( *this ) , detail::_1 ) );
  111. sys( in , m_dxdt.m_v , t );
  112. detail::adams_moulton_call_algebra< steps , algebra_type , operations_type >()( m_algebra , in , out , m_dxdt.m_v , buf , m_coefficients , dt );
  113. }
  114. template< class StateType >
  115. void adjust_size( const StateType &x )
  116. {
  117. resize_impl( x );
  118. }
  119. algebra_type& algebra()
  120. { return m_algebra; }
  121. const algebra_type& algebra() const
  122. { return m_algebra; }
  123. private:
  124. template< class StateIn >
  125. bool resize_impl( const StateIn &x )
  126. {
  127. return adjust_size_by_resizeability( m_dxdt , x , typename is_resizeable<deriv_type>::type() );
  128. }
  129. const detail::adams_moulton_coefficients< value_type , steps > m_coefficients;
  130. wrapped_deriv_type m_dxdt;
  131. resizer_type m_resizer;
  132. protected:
  133. algebra_type m_algebra_instance;
  134. algebra_type &m_algebra;
  135. };
  136. } // odeint
  137. } // numeric
  138. } // boost
  139. #endif // BOOST_NUMERIC_ODEINT_STEPPER_ADAMS_MOULTON_HPP_INCLUDED