stepper_details.cpp 5.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200
  1. /*
  2. * stepper_details.cpp
  3. *
  4. * This example demonstrates some details about the steppers in odeint.
  5. *
  6. *
  7. * Copyright 2011-2012 Karsten Ahnert
  8. * Copyright 2012 Mario Mulansky
  9. * Copyright 2013 Pascal Germroth
  10. *
  11. * Distributed under the Boost Software License, Version 1.0.
  12. * (See accompanying file LICENSE_1_0.txt or
  13. * copy at http://www.boost.org/LICENSE_1_0.txt)
  14. */
  15. #include <iostream>
  16. #include <boost/array.hpp>
  17. #include <boost/bind.hpp>
  18. #include <boost/numeric/odeint.hpp>
  19. using namespace std;
  20. using namespace boost::numeric::odeint;
  21. const size_t N = 3;
  22. typedef boost::array< double , N > state_type;
  23. //[ system_function_structure
  24. void sys( const state_type & /*x*/ , state_type & /*dxdt*/ , const double /*t*/ )
  25. {
  26. // ...
  27. }
  28. //]
  29. void sys1( const state_type &/*x*/ , state_type &/*dxdt*/ , const double /*t*/ )
  30. {
  31. }
  32. void sys2( const state_type &/*x*/ , state_type &/*dxdt*/ , const double /*t*/ )
  33. {
  34. }
  35. //[ symplectic_stepper_detail_system_function
  36. typedef boost::array< double , 1 > vector_type;
  37. struct harm_osc_f1
  38. {
  39. void operator()( const vector_type &p , vector_type &dqdt )
  40. {
  41. dqdt[0] = p[0];
  42. }
  43. };
  44. struct harm_osc_f2
  45. {
  46. void operator()( const vector_type &q , vector_type &dpdt )
  47. {
  48. dpdt[0] = -q[0];
  49. }
  50. };
  51. //]
  52. //[ symplectic_stepper_detail_system_class
  53. struct harm_osc
  54. {
  55. void f1( const vector_type &p , vector_type &dqdt ) const
  56. {
  57. dqdt[0] = p[0];
  58. }
  59. void f2( const vector_type &q , vector_type &dpdt ) const
  60. {
  61. dpdt[0] = -q[0];
  62. }
  63. };
  64. //]
  65. int main( int argc , char **argv )
  66. {
  67. using namespace std;
  68. // Explicit stepper example
  69. {
  70. double t( 0.0 ) , dt( 0.1 );
  71. state_type in , out , dxdtin , inout;
  72. //[ explicit_stepper_detail_example
  73. runge_kutta4< state_type > rk;
  74. rk.do_step( sys1 , inout , t , dt ); // In-place transformation of inout
  75. rk.do_step( sys2 , inout , t , dt ); // call with different system: Ok
  76. rk.do_step( sys1 , in , t , out , dt ); // Out-of-place transformation
  77. rk.do_step( sys1 , inout , dxdtin , t , dt ); // In-place tranformation of inout
  78. rk.do_step( sys1 , in , dxdtin , t , out , dt ); // Out-of-place transformation
  79. //]
  80. }
  81. // FSAL stepper example
  82. {
  83. double t( 0.0 ) , dt( 0.1 );
  84. state_type in , in2 , in3 , out , dxdtin , dxdtout , inout , dxdtinout;
  85. //[ fsal_stepper_detail_example
  86. runge_kutta_dopri5< state_type > rk;
  87. rk.do_step( sys1 , in , t , out , dt );
  88. rk.do_step( sys2 , in , t , out , dt ); // DONT do this, sys1 is assumed
  89. rk.do_step( sys2 , in2 , t , out , dt );
  90. rk.do_step( sys2 , in3 , t , out , dt ); // DONT do this, in2 is assumed
  91. rk.do_step( sys1 , inout , dxdtinout , t , dt );
  92. rk.do_step( sys2 , inout , dxdtinout , t , dt ); // Ok, internal derivative is not used, dxdtinout is updated
  93. rk.do_step( sys1 , in , dxdtin , t , out , dxdtout , dt );
  94. rk.do_step( sys2 , in , dxdtin , t , out , dxdtout , dt ); // Ok, internal derivative is not used
  95. //]
  96. }
  97. // Symplectic harmonic oscillator example
  98. {
  99. double t( 0.0 ) , dt( 0.1 );
  100. //[ symplectic_stepper_detail_example
  101. pair< vector_type , vector_type > x;
  102. x.first[0] = 1.0; x.second[0] = 0.0;
  103. symplectic_rkn_sb3a_mclachlan< vector_type > rkn;
  104. rkn.do_step( make_pair( harm_osc_f1() , harm_osc_f2() ) , x , t , dt );
  105. //]
  106. //[ symplectic_stepper_detail_system_class_example
  107. harm_osc h;
  108. rkn.do_step( make_pair( boost::bind( &harm_osc::f1 , h , _1 , _2 ) , boost::bind( &harm_osc::f2 , h , _1 , _2 ) ) ,
  109. x , t , dt );
  110. //]
  111. }
  112. // Simplified harmonic oscillator example
  113. {
  114. double t = 0.0, dt = 0.1;
  115. //[ simplified_symplectic_stepper_example
  116. pair< vector_type , vector_type > x;
  117. x.first[0] = 1.0; x.second[0] = 0.0;
  118. symplectic_rkn_sb3a_mclachlan< vector_type > rkn;
  119. rkn.do_step( harm_osc_f1() , x , t , dt );
  120. //]
  121. vector_type q = {{ 1.0 }} , p = {{ 0.0 }};
  122. //[ symplectic_stepper_detail_ref_usage
  123. rkn.do_step( harm_osc_f1() , make_pair( boost::ref( q ) , boost::ref( p ) ) , t , dt );
  124. rkn.do_step( harm_osc_f1() , q , p , t , dt );
  125. rkn.do_step( make_pair( harm_osc_f1() , harm_osc_f2() ) , q , p , t , dt );
  126. //]
  127. }
  128. // adams_bashforth_moulton stepper example
  129. {
  130. double t = 0.0 , dt = 0.1;
  131. state_type inout;
  132. //[ multistep_detail_example
  133. adams_bashforth_moulton< 5 , state_type > abm;
  134. abm.initialize( sys , inout , t , dt );
  135. abm.do_step( sys , inout , t , dt );
  136. //]
  137. //[ multistep_detail_own_stepper_initialization
  138. abm.initialize( runge_kutta_fehlberg78< state_type >() , sys , inout , t , dt );
  139. //]
  140. }
  141. // dense output stepper examples
  142. {
  143. double t = 0.0 , dt = 0.1;
  144. state_type in;
  145. //[ dense_output_detail_example
  146. dense_output_runge_kutta< controlled_runge_kutta< runge_kutta_dopri5< state_type > > > dense;
  147. dense.initialize( in , t , dt );
  148. pair< double , double > times = dense.do_step( sys );
  149. (void)times;
  150. //]
  151. state_type inout;
  152. double t_start = 0.0 , t_end = 1.0;
  153. //[ dense_output_detail_generation1
  154. typedef boost::numeric::odeint::result_of::make_dense_output<
  155. runge_kutta_dopri5< state_type > >::type dense_stepper_type;
  156. dense_stepper_type dense2 = make_dense_output( 1.0e-6 , 1.0e-6 , runge_kutta_dopri5< state_type >() );
  157. (void)dense2;
  158. //]
  159. //[ dense_output_detail_generation2
  160. integrate_const( make_dense_output( 1.0e-6 , 1.0e-6 , runge_kutta_dopri5< state_type >() ) , sys , inout , t_start , t_end , dt );
  161. //]
  162. }
  163. return 0;
  164. }