rosenbrock4.hpp 14 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345
  1. /*
  2. [auto_generated]
  3. boost/numeric/odeint/stepper/rosenbrock4.hpp
  4. [begin_description]
  5. Implementation of the Rosenbrock 4 method for solving stiff ODEs. Note, that a
  6. controller and a dense-output stepper exist for this 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_ROSENBROCK4_HPP_INCLUDED
  15. #define BOOST_NUMERIC_ODEINT_STEPPER_ROSENBROCK4_HPP_INCLUDED
  16. #include <boost/numeric/odeint/util/bind.hpp>
  17. #include <boost/numeric/odeint/util/unwrap_reference.hpp>
  18. #include <boost/numeric/ublas/vector.hpp>
  19. #include <boost/numeric/ublas/matrix.hpp>
  20. #include <boost/numeric/ublas/lu.hpp>
  21. #include <boost/numeric/odeint/stepper/stepper_categories.hpp>
  22. #include <boost/numeric/odeint/util/ublas_wrapper.hpp>
  23. #include <boost/numeric/odeint/util/is_resizeable.hpp>
  24. #include <boost/numeric/odeint/util/resizer.hpp>
  25. #include <boost/numeric/ublas/vector.hpp>
  26. #include <boost/numeric/ublas/matrix.hpp>
  27. #include <boost/numeric/ublas/lu.hpp>
  28. namespace boost {
  29. namespace numeric {
  30. namespace odeint {
  31. /*
  32. * ToDo:
  33. *
  34. * 2. Interfacing for odeint, check if controlled_error_stepper can be used
  35. * 3. dense output
  36. */
  37. template< class Value >
  38. struct default_rosenbrock_coefficients
  39. {
  40. typedef Value value_type;
  41. typedef unsigned short order_type;
  42. default_rosenbrock_coefficients( void )
  43. : gamma ( static_cast< value_type >( 0.25 ) ) ,
  44. d1 ( static_cast< value_type >( 0.25 ) ) ,
  45. d2 ( static_cast< value_type >( -0.1043 ) ) ,
  46. d3 ( static_cast< value_type >( 0.1035 ) ) ,
  47. d4 ( static_cast< value_type >( 0.3620000000000023e-01 ) ) ,
  48. c2 ( static_cast< value_type >( 0.386 ) ) ,
  49. c3 ( static_cast< value_type >( 0.21 ) ) ,
  50. c4 ( static_cast< value_type >( 0.63 ) ) ,
  51. c21 ( static_cast< value_type >( -0.5668800000000000e+01 ) ) ,
  52. a21 ( static_cast< value_type >( 0.1544000000000000e+01 ) ) ,
  53. c31 ( static_cast< value_type >( -0.2430093356833875e+01 ) ) ,
  54. c32 ( static_cast< value_type >( -0.2063599157091915e+00 ) ) ,
  55. a31 ( static_cast< value_type >( 0.9466785280815826e+00 ) ) ,
  56. a32 ( static_cast< value_type >( 0.2557011698983284e+00 ) ) ,
  57. c41 ( static_cast< value_type >( -0.1073529058151375e+00 ) ) ,
  58. c42 ( static_cast< value_type >( -0.9594562251023355e+01 ) ) ,
  59. c43 ( static_cast< value_type >( -0.2047028614809616e+02 ) ) ,
  60. a41 ( static_cast< value_type >( 0.3314825187068521e+01 ) ) ,
  61. a42 ( static_cast< value_type >( 0.2896124015972201e+01 ) ) ,
  62. a43 ( static_cast< value_type >( 0.9986419139977817e+00 ) ) ,
  63. c51 ( static_cast< value_type >( 0.7496443313967647e+01 ) ) ,
  64. c52 ( static_cast< value_type >( -0.1024680431464352e+02 ) ) ,
  65. c53 ( static_cast< value_type >( -0.3399990352819905e+02 ) ) ,
  66. c54 ( static_cast< value_type >( 0.1170890893206160e+02 ) ) ,
  67. a51 ( static_cast< value_type >( 0.1221224509226641e+01 ) ) ,
  68. a52 ( static_cast< value_type >( 0.6019134481288629e+01 ) ) ,
  69. a53 ( static_cast< value_type >( 0.1253708332932087e+02 ) ) ,
  70. a54 ( static_cast< value_type >( -0.6878860361058950e+00 ) ) ,
  71. c61 ( static_cast< value_type >( 0.8083246795921522e+01 ) ) ,
  72. c62 ( static_cast< value_type >( -0.7981132988064893e+01 ) ) ,
  73. c63 ( static_cast< value_type >( -0.3152159432874371e+02 ) ) ,
  74. c64 ( static_cast< value_type >( 0.1631930543123136e+02 ) ) ,
  75. c65 ( static_cast< value_type >( -0.6058818238834054e+01 ) ) ,
  76. d21 ( static_cast< value_type >( 0.1012623508344586e+02 ) ) ,
  77. d22 ( static_cast< value_type >( -0.7487995877610167e+01 ) ) ,
  78. d23 ( static_cast< value_type >( -0.3480091861555747e+02 ) ) ,
  79. d24 ( static_cast< value_type >( -0.7992771707568823e+01 ) ) ,
  80. d25 ( static_cast< value_type >( 0.1025137723295662e+01 ) ) ,
  81. d31 ( static_cast< value_type >( -0.6762803392801253e+00 ) ) ,
  82. d32 ( static_cast< value_type >( 0.6087714651680015e+01 ) ) ,
  83. d33 ( static_cast< value_type >( 0.1643084320892478e+02 ) ) ,
  84. d34 ( static_cast< value_type >( 0.2476722511418386e+02 ) ) ,
  85. d35 ( static_cast< value_type >( -0.6594389125716872e+01 ) )
  86. {}
  87. const value_type gamma;
  88. const value_type d1 , d2 , d3 , d4;
  89. const value_type c2 , c3 , c4;
  90. const value_type c21 ;
  91. const value_type a21;
  92. const value_type c31 , c32;
  93. const value_type a31 , a32;
  94. const value_type c41 , c42 , c43;
  95. const value_type a41 , a42 , a43;
  96. const value_type c51 , c52 , c53 , c54;
  97. const value_type a51 , a52 , a53 , a54;
  98. const value_type c61 , c62 , c63 , c64 , c65;
  99. const value_type d21 , d22 , d23 , d24 , d25;
  100. const value_type d31 , d32 , d33 , d34 , d35;
  101. static const order_type stepper_order = 4;
  102. static const order_type error_order = 3;
  103. };
  104. template< class Value , class Coefficients = default_rosenbrock_coefficients< Value > , class Resizer = initially_resizer >
  105. class rosenbrock4
  106. {
  107. private:
  108. public:
  109. typedef Value value_type;
  110. typedef boost::numeric::ublas::vector< value_type > state_type;
  111. typedef state_type deriv_type;
  112. typedef value_type time_type;
  113. typedef boost::numeric::ublas::matrix< value_type > matrix_type;
  114. typedef boost::numeric::ublas::permutation_matrix< size_t > pmatrix_type;
  115. typedef Resizer resizer_type;
  116. typedef Coefficients rosenbrock_coefficients;
  117. typedef stepper_tag stepper_category;
  118. typedef unsigned short order_type;
  119. typedef state_wrapper< state_type > wrapped_state_type;
  120. typedef state_wrapper< deriv_type > wrapped_deriv_type;
  121. typedef state_wrapper< matrix_type > wrapped_matrix_type;
  122. typedef state_wrapper< pmatrix_type > wrapped_pmatrix_type;
  123. typedef rosenbrock4< Value , Coefficients , Resizer > stepper_type;
  124. const static order_type stepper_order = rosenbrock_coefficients::stepper_order;
  125. const static order_type error_order = rosenbrock_coefficients::error_order;
  126. rosenbrock4( void )
  127. : m_resizer() , m_x_err_resizer() ,
  128. m_jac() , m_pm() ,
  129. m_dfdt() , m_dxdt() , m_dxdtnew() ,
  130. m_g1() , m_g2() , m_g3() , m_g4() , m_g5() ,
  131. m_cont3() , m_cont4() , m_xtmp() , m_x_err() ,
  132. m_coef()
  133. { }
  134. order_type order() const { return stepper_order; }
  135. template< class System >
  136. void do_step( System system , const state_type &x , time_type t , state_type &xout , time_type dt , state_type &xerr )
  137. {
  138. // get the system and jacobi function
  139. typedef typename odeint::unwrap_reference< System >::type system_type;
  140. typedef typename odeint::unwrap_reference< typename system_type::first_type >::type deriv_func_type;
  141. typedef typename odeint::unwrap_reference< typename system_type::second_type >::type jacobi_func_type;
  142. system_type &sys = system;
  143. deriv_func_type &deriv_func = sys.first;
  144. jacobi_func_type &jacobi_func = sys.second;
  145. const size_t n = x.size();
  146. m_resizer.adjust_size( x , detail::bind( &stepper_type::template resize_impl<state_type> , detail::ref( *this ) , detail::_1 ) );
  147. for( size_t i=0 ; i<n ; ++i )
  148. m_pm.m_v( i ) = i;
  149. deriv_func( x , m_dxdt.m_v , t );
  150. jacobi_func( x , m_jac.m_v , t , m_dfdt.m_v );
  151. m_jac.m_v *= -1.0;
  152. m_jac.m_v += 1.0 / m_coef.gamma / dt * boost::numeric::ublas::identity_matrix< value_type >( n );
  153. boost::numeric::ublas::lu_factorize( m_jac.m_v , m_pm.m_v );
  154. for( size_t i=0 ; i<n ; ++i )
  155. m_g1.m_v[i] = m_dxdt.m_v[i] + dt * m_coef.d1 * m_dfdt.m_v[i];
  156. boost::numeric::ublas::lu_substitute( m_jac.m_v , m_pm.m_v , m_g1.m_v );
  157. for( size_t i=0 ; i<n ; ++i )
  158. m_xtmp.m_v[i] = x[i] + m_coef.a21 * m_g1.m_v[i];
  159. deriv_func( m_xtmp.m_v , m_dxdtnew.m_v , t + m_coef.c2 * dt );
  160. for( size_t i=0 ; i<n ; ++i )
  161. m_g2.m_v[i] = m_dxdtnew.m_v[i] + dt * m_coef.d2 * m_dfdt.m_v[i] + m_coef.c21 * m_g1.m_v[i] / dt;
  162. boost::numeric::ublas::lu_substitute( m_jac.m_v , m_pm.m_v , m_g2.m_v );
  163. for( size_t i=0 ; i<n ; ++i )
  164. m_xtmp.m_v[i] = x[i] + m_coef.a31 * m_g1.m_v[i] + m_coef.a32 * m_g2.m_v[i];
  165. deriv_func( m_xtmp.m_v , m_dxdtnew.m_v , t + m_coef.c3 * dt );
  166. for( size_t i=0 ; i<n ; ++i )
  167. m_g3.m_v[i] = m_dxdtnew.m_v[i] + dt * m_coef.d3 * m_dfdt.m_v[i] + ( m_coef.c31 * m_g1.m_v[i] + m_coef.c32 * m_g2.m_v[i] ) / dt;
  168. boost::numeric::ublas::lu_substitute( m_jac.m_v , m_pm.m_v , m_g3.m_v );
  169. for( size_t i=0 ; i<n ; ++i )
  170. m_xtmp.m_v[i] = x[i] + m_coef.a41 * m_g1.m_v[i] + m_coef.a42 * m_g2.m_v[i] + m_coef.a43 * m_g3.m_v[i];
  171. deriv_func( m_xtmp.m_v , m_dxdtnew.m_v , t + m_coef.c4 * dt );
  172. for( size_t i=0 ; i<n ; ++i )
  173. m_g4.m_v[i] = m_dxdtnew.m_v[i] + dt * m_coef.d4 * m_dfdt.m_v[i] + ( m_coef.c41 * m_g1.m_v[i] + m_coef.c42 * m_g2.m_v[i] + m_coef.c43 * m_g3.m_v[i] ) / dt;
  174. boost::numeric::ublas::lu_substitute( m_jac.m_v , m_pm.m_v , m_g4.m_v );
  175. for( size_t i=0 ; i<n ; ++i )
  176. m_xtmp.m_v[i] = x[i] + m_coef.a51 * m_g1.m_v[i] + m_coef.a52 * m_g2.m_v[i] + m_coef.a53 * m_g3.m_v[i] + m_coef.a54 * m_g4.m_v[i];
  177. deriv_func( m_xtmp.m_v , m_dxdtnew.m_v , t + dt );
  178. for( size_t i=0 ; i<n ; ++i )
  179. m_g5.m_v[i] = m_dxdtnew.m_v[i] + ( m_coef.c51 * m_g1.m_v[i] + m_coef.c52 * m_g2.m_v[i] + m_coef.c53 * m_g3.m_v[i] + m_coef.c54 * m_g4.m_v[i] ) / dt;
  180. boost::numeric::ublas::lu_substitute( m_jac.m_v , m_pm.m_v , m_g5.m_v );
  181. for( size_t i=0 ; i<n ; ++i )
  182. m_xtmp.m_v[i] += m_g5.m_v[i];
  183. deriv_func( m_xtmp.m_v , m_dxdtnew.m_v , t + dt );
  184. for( size_t i=0 ; i<n ; ++i )
  185. xerr[i] = m_dxdtnew.m_v[i] + ( m_coef.c61 * m_g1.m_v[i] + m_coef.c62 * m_g2.m_v[i] + m_coef.c63 * m_g3.m_v[i] + m_coef.c64 * m_g4.m_v[i] + m_coef.c65 * m_g5.m_v[i] ) / dt;
  186. boost::numeric::ublas::lu_substitute( m_jac.m_v , m_pm.m_v , xerr );
  187. for( size_t i=0 ; i<n ; ++i )
  188. xout[i] = m_xtmp.m_v[i] + xerr[i];
  189. }
  190. template< class System >
  191. void do_step( System system , state_type &x , time_type t , time_type dt , state_type &xerr )
  192. {
  193. do_step( system , x , t , x , dt , xerr );
  194. }
  195. /*
  196. * do_step without error output - just calls above functions with and neglects the error estimate
  197. */
  198. template< class System >
  199. void do_step( System system , const state_type &x , time_type t , state_type &xout , time_type dt )
  200. {
  201. m_x_err_resizer.adjust_size( x , detail::bind( &stepper_type::template resize_x_err<state_type> , detail::ref( *this ) , detail::_1 ) );
  202. do_step( system , x , t , xout , dt , m_x_err.m_v );
  203. }
  204. template< class System >
  205. void do_step( System system , state_type &x , time_type t , time_type dt )
  206. {
  207. m_x_err_resizer.adjust_size( x , detail::bind( &stepper_type::template resize_x_err<state_type> , detail::ref( *this ) , detail::_1 ) );
  208. do_step( system , x , t , dt , m_x_err.m_v );
  209. }
  210. void prepare_dense_output()
  211. {
  212. const size_t n = m_g1.m_v.size();
  213. for( size_t i=0 ; i<n ; ++i )
  214. {
  215. m_cont3.m_v[i] = m_coef.d21 * m_g1.m_v[i] + m_coef.d22 * m_g2.m_v[i] + m_coef.d23 * m_g3.m_v[i] + m_coef.d24 * m_g4.m_v[i] + m_coef.d25 * m_g5.m_v[i];
  216. m_cont4.m_v[i] = m_coef.d31 * m_g1.m_v[i] + m_coef.d32 * m_g2.m_v[i] + m_coef.d33 * m_g3.m_v[i] + m_coef.d34 * m_g4.m_v[i] + m_coef.d35 * m_g5.m_v[i];
  217. }
  218. }
  219. void calc_state( time_type t , state_type &x ,
  220. const state_type &x_old , time_type t_old ,
  221. const state_type &x_new , time_type t_new )
  222. {
  223. const size_t n = m_g1.m_v.size();
  224. time_type dt = t_new - t_old;
  225. time_type s = ( t - t_old ) / dt;
  226. time_type s1 = 1.0 - s;
  227. for( size_t i=0 ; i<n ; ++i )
  228. x[i] = x_old[i] * s1 + s * ( x_new[i] + s1 * ( m_cont3.m_v[i] + s * m_cont4.m_v[i] ) );
  229. }
  230. template< class StateType >
  231. void adjust_size( const StateType &x )
  232. {
  233. resize_impl( x );
  234. resize_x_err( x );
  235. }
  236. protected:
  237. template< class StateIn >
  238. bool resize_impl( const StateIn &x )
  239. {
  240. bool resized = false;
  241. resized |= adjust_size_by_resizeability( m_dxdt , x , typename is_resizeable<deriv_type>::type() );
  242. resized |= adjust_size_by_resizeability( m_dfdt , x , typename is_resizeable<deriv_type>::type() );
  243. resized |= adjust_size_by_resizeability( m_dxdtnew , x , typename is_resizeable<deriv_type>::type() );
  244. resized |= adjust_size_by_resizeability( m_xtmp , x , typename is_resizeable<state_type>::type() );
  245. resized |= adjust_size_by_resizeability( m_g1 , x , typename is_resizeable<state_type>::type() );
  246. resized |= adjust_size_by_resizeability( m_g2 , x , typename is_resizeable<state_type>::type() );
  247. resized |= adjust_size_by_resizeability( m_g3 , x , typename is_resizeable<state_type>::type() );
  248. resized |= adjust_size_by_resizeability( m_g4 , x , typename is_resizeable<state_type>::type() );
  249. resized |= adjust_size_by_resizeability( m_g5 , x , typename is_resizeable<state_type>::type() );
  250. resized |= adjust_size_by_resizeability( m_cont3 , x , typename is_resizeable<state_type>::type() );
  251. resized |= adjust_size_by_resizeability( m_cont4 , x , typename is_resizeable<state_type>::type() );
  252. resized |= adjust_size_by_resizeability( m_jac , x , typename is_resizeable<matrix_type>::type() );
  253. resized |= adjust_size_by_resizeability( m_pm , x , typename is_resizeable<pmatrix_type>::type() );
  254. return resized;
  255. }
  256. template< class StateIn >
  257. bool resize_x_err( const StateIn &x )
  258. {
  259. return adjust_size_by_resizeability( m_x_err , x , typename is_resizeable<state_type>::type() );
  260. }
  261. private:
  262. resizer_type m_resizer;
  263. resizer_type m_x_err_resizer;
  264. wrapped_matrix_type m_jac;
  265. wrapped_pmatrix_type m_pm;
  266. wrapped_deriv_type m_dfdt , m_dxdt , m_dxdtnew;
  267. wrapped_state_type m_g1 , m_g2 , m_g3 , m_g4 , m_g5;
  268. wrapped_state_type m_cont3 , m_cont4;
  269. wrapped_state_type m_xtmp;
  270. wrapped_state_type m_x_err;
  271. const rosenbrock_coefficients m_coef;
  272. };
  273. } // namespace odeint
  274. } // namespace numeric
  275. } // namespace boost
  276. #endif // BOOST_NUMERIC_ODEINT_STEPPER_ROSENBROCK4_HPP_INCLUDED