generic_rk_algorithm.hpp 8.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247
  1. /*
  2. [auto_generated]
  3. boost/numeric/odeint/stepper/detail/generic_rk_algorithm.hpp
  4. [begin_description]
  5. Implementation of the generic Runge-Kutta method.
  6. [end_description]
  7. Copyright 2011-2013 Mario Mulansky
  8. Copyright 2011-2012 Karsten Ahnert
  9. Copyright 2012 Christoph Koke
  10. Distributed under the Boost Software License, Version 1.0.
  11. (See accompanying file LICENSE_1_0.txt or
  12. copy at http://www.boost.org/LICENSE_1_0.txt)
  13. */
  14. #ifndef BOOST_NUMERIC_ODEINT_STEPPER_DETAIL_GENERIC_RK_ALGORITHM_HPP_INCLUDED
  15. #define BOOST_NUMERIC_ODEINT_STEPPER_DETAIL_GENERIC_RK_ALGORITHM_HPP_INCLUDED
  16. #include <boost/static_assert.hpp>
  17. #include <boost/mpl/vector.hpp>
  18. #include <boost/mpl/push_back.hpp>
  19. #include <boost/mpl/for_each.hpp>
  20. #include <boost/mpl/range_c.hpp>
  21. #include <boost/mpl/copy.hpp>
  22. #include <boost/mpl/size_t.hpp>
  23. #include <boost/fusion/algorithm.hpp>
  24. #include <boost/fusion/iterator.hpp>
  25. #include <boost/fusion/mpl.hpp>
  26. #include <boost/fusion/sequence.hpp>
  27. #include <boost/array.hpp>
  28. #include <boost/numeric/odeint/algebra/range_algebra.hpp>
  29. #include <boost/numeric/odeint/algebra/default_operations.hpp>
  30. #include <boost/numeric/odeint/stepper/detail/generic_rk_call_algebra.hpp>
  31. #include <boost/numeric/odeint/stepper/detail/generic_rk_operations.hpp>
  32. #include <boost/numeric/odeint/util/bind.hpp>
  33. namespace boost {
  34. namespace numeric {
  35. namespace odeint {
  36. namespace detail {
  37. template< class T , class Constant >
  38. struct array_wrapper
  39. {
  40. typedef const typename boost::array< T , Constant::value > type;
  41. };
  42. template< class T , size_t i >
  43. struct stage
  44. {
  45. T c;
  46. boost::array< T , i > a;
  47. };
  48. template< class T , class Constant >
  49. struct stage_wrapper
  50. {
  51. typedef stage< T , Constant::value > type;
  52. };
  53. template<
  54. size_t StageCount,
  55. class Value ,
  56. class Algebra ,
  57. class Operations
  58. >
  59. class generic_rk_algorithm {
  60. public:
  61. typedef mpl::range_c< size_t , 1 , StageCount > stage_indices;
  62. typedef typename boost::fusion::result_of::as_vector
  63. <
  64. typename boost::mpl::copy
  65. <
  66. stage_indices ,
  67. boost::mpl::inserter
  68. <
  69. boost::mpl::vector0< > ,
  70. boost::mpl::push_back< boost::mpl::_1 , array_wrapper< Value , boost::mpl::_2 > >
  71. >
  72. >::type
  73. >::type coef_a_type;
  74. typedef boost::array< Value , StageCount > coef_b_type;
  75. typedef boost::array< Value , StageCount > coef_c_type;
  76. typedef typename boost::fusion::result_of::as_vector
  77. <
  78. typename boost::mpl::push_back
  79. <
  80. typename boost::mpl::copy
  81. <
  82. stage_indices,
  83. boost::mpl::inserter
  84. <
  85. boost::mpl::vector0<> ,
  86. boost::mpl::push_back< boost::mpl::_1 , stage_wrapper< Value , boost::mpl::_2 > >
  87. >
  88. >::type ,
  89. stage< Value , StageCount >
  90. >::type
  91. >::type stage_vector_base;
  92. struct stage_vector : public stage_vector_base
  93. {
  94. struct do_insertion
  95. {
  96. stage_vector_base &m_base;
  97. const coef_a_type &m_a;
  98. const coef_c_type &m_c;
  99. do_insertion( stage_vector_base &base , const coef_a_type &a , const coef_c_type &c )
  100. : m_base( base ) , m_a( a ) , m_c( c ) { }
  101. template< class Index >
  102. void operator()( Index ) const
  103. {
  104. //boost::fusion::at< Index >( m_base ) = stage< double , Index::value+1 , intermediate_stage >( m_c[ Index::value ] , boost::fusion::at< Index >( m_a ) );
  105. boost::fusion::at< Index >( m_base ).c = m_c[ Index::value ];
  106. boost::fusion::at< Index >( m_base ).a = boost::fusion::at< Index >( m_a );
  107. }
  108. };
  109. struct print_butcher
  110. {
  111. const stage_vector_base &m_base;
  112. std::ostream &m_os;
  113. print_butcher( const stage_vector_base &base , std::ostream &os )
  114. : m_base( base ) , m_os( os )
  115. { }
  116. template<class Index>
  117. void operator()(Index) const {
  118. m_os << boost::fusion::at<Index>(m_base).c << " | ";
  119. for( size_t i=0 ; i<Index::value ; ++i )
  120. m_os << boost::fusion::at<Index>(m_base).a[i] << " ";
  121. m_os << std::endl;
  122. }
  123. };
  124. stage_vector( const coef_a_type &a , const coef_b_type &b , const coef_c_type &c )
  125. {
  126. typedef boost::mpl::range_c< size_t , 0 , StageCount-1 > indices;
  127. boost::mpl::for_each< indices >( do_insertion( *this , a , c ) );
  128. boost::fusion::at_c< StageCount - 1 >( *this ).c = c[ StageCount - 1 ];
  129. boost::fusion::at_c< StageCount - 1 >( *this ).a = b;
  130. }
  131. void print( std::ostream &os ) const
  132. {
  133. typedef boost::mpl::range_c< size_t , 0 , StageCount > indices;
  134. boost::mpl::for_each< indices >( print_butcher( *this , os ) );
  135. }
  136. };
  137. template< class System , class StateIn , class StateTemp , class DerivIn , class Deriv ,
  138. class StateOut , class Time >
  139. struct calculate_stage
  140. {
  141. Algebra &algebra;
  142. System &system;
  143. const StateIn &x;
  144. StateTemp &x_tmp;
  145. StateOut &x_out;
  146. const DerivIn &dxdt;
  147. Deriv *F;
  148. Time t;
  149. Time dt;
  150. calculate_stage( Algebra &_algebra , System &_system , const StateIn &_x , const DerivIn &_dxdt , StateOut &_out ,
  151. StateTemp &_x_tmp , Deriv *_F , Time _t , Time _dt )
  152. : algebra( _algebra ) , system( _system ) , x( _x ) , x_tmp( _x_tmp ) , x_out( _out) , dxdt( _dxdt ) , F( _F ) , t( _t ) , dt( _dt )
  153. {}
  154. template< typename T , size_t stage_number >
  155. void inline operator()( stage< T , stage_number > const &stage ) const
  156. //typename stage_fusion_wrapper< T , mpl::size_t< stage_number > , intermediate_stage >::type const &stage ) const
  157. {
  158. if( stage_number > 1 )
  159. {
  160. #ifdef BOOST_MSVC
  161. #pragma warning( disable : 4307 34 )
  162. #endif
  163. system( x_tmp , F[stage_number-2].m_v , t + stage.c * dt );
  164. #ifdef BOOST_MSVC
  165. #pragma warning( default : 4307 34 )
  166. #endif
  167. }
  168. //std::cout << stage_number-2 << ", t': " << t + stage.c * dt << std::endl;
  169. if( stage_number < StageCount )
  170. detail::template generic_rk_call_algebra< stage_number , Algebra >()( algebra , x_tmp , x , dxdt , F ,
  171. detail::generic_rk_scale_sum< stage_number , Operations , Value , Time >( stage.a , dt) );
  172. // algebra_type::template for_eachn<stage_number>( x_tmp , x , dxdt , F ,
  173. // typename operations_type::template scale_sumn< stage_number , time_type >( stage.a , dt ) );
  174. else
  175. detail::template generic_rk_call_algebra< stage_number , Algebra >()( algebra , x_out , x , dxdt , F ,
  176. detail::generic_rk_scale_sum< stage_number , Operations , Value , Time >( stage.a , dt ) );
  177. // algebra_type::template for_eachn<stage_number>( x_out , x , dxdt , F ,
  178. // typename operations_type::template scale_sumn< stage_number , time_type >( stage.a , dt ) );
  179. }
  180. };
  181. generic_rk_algorithm( const coef_a_type &a , const coef_b_type &b , const coef_c_type &c )
  182. : m_stages( a , b , c )
  183. { }
  184. template< class System , class StateIn , class DerivIn , class Time , class StateOut , class StateTemp , class Deriv >
  185. void inline do_step( Algebra &algebra , System system , const StateIn &in , const DerivIn &dxdt ,
  186. Time t , StateOut &out , Time dt ,
  187. StateTemp &x_tmp , Deriv F[StageCount-1] ) const
  188. {
  189. typedef typename odeint::unwrap_reference< System >::type unwrapped_system_type;
  190. unwrapped_system_type &sys = system;
  191. boost::fusion::for_each( m_stages , calculate_stage<
  192. unwrapped_system_type , StateIn , StateTemp , DerivIn , Deriv , StateOut , Time >
  193. ( algebra , sys , in , dxdt , out , x_tmp , F , t , dt ) );
  194. }
  195. private:
  196. stage_vector m_stages;
  197. };
  198. }
  199. }
  200. }
  201. }
  202. #endif // BOOST_NUMERIC_ODEINT_STEPPER_DETAIL_GENERIC_RK_ALGORITHM_HPP_INCLUDED