harmonic_oscillator.cpp 4.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211
  1. /*
  2. Copyright 2010-2012 Karsten Ahnert
  3. Copyright 2011-2013 Mario Mulansky
  4. Copyright 2013 Pascal Germroth
  5. Distributed under the Boost Software License, Version 1.0.
  6. (See accompanying file LICENSE_1_0.txt or
  7. copy at http://www.boost.org/LICENSE_1_0.txt)
  8. */
  9. #include <iostream>
  10. #include <vector>
  11. #include <boost/numeric/odeint.hpp>
  12. //[ rhs_function
  13. /* The type of container used to hold the state vector */
  14. typedef std::vector< double > state_type;
  15. const double gam = 0.15;
  16. /* The rhs of x' = f(x) */
  17. void harmonic_oscillator( const state_type &x , state_type &dxdt , const double /* t */ )
  18. {
  19. dxdt[0] = x[1];
  20. dxdt[1] = -x[0] - gam*x[1];
  21. }
  22. //]
  23. //[ rhs_class
  24. /* The rhs of x' = f(x) defined as a class */
  25. class harm_osc {
  26. double m_gam;
  27. public:
  28. harm_osc( double gam ) : m_gam(gam) { }
  29. void operator() ( const state_type &x , state_type &dxdt , const double /* t */ )
  30. {
  31. dxdt[0] = x[1];
  32. dxdt[1] = -x[0] - m_gam*x[1];
  33. }
  34. };
  35. //]
  36. //[ integrate_observer
  37. struct push_back_state_and_time
  38. {
  39. std::vector< state_type >& m_states;
  40. std::vector< double >& m_times;
  41. push_back_state_and_time( std::vector< state_type > &states , std::vector< double > &times )
  42. : m_states( states ) , m_times( times ) { }
  43. void operator()( const state_type &x , double t )
  44. {
  45. m_states.push_back( x );
  46. m_times.push_back( t );
  47. }
  48. };
  49. //]
  50. struct write_state
  51. {
  52. void operator()( const state_type &x ) const
  53. {
  54. std::cout << x[0] << "\t" << x[1] << "\n";
  55. }
  56. };
  57. int main(int /* argc */ , char** /* argv */ )
  58. {
  59. using namespace std;
  60. using namespace boost::numeric::odeint;
  61. //[ state_initialization
  62. state_type x(2);
  63. x[0] = 1.0; // start at x=1.0, p=0.0
  64. x[1] = 0.0;
  65. //]
  66. //[ integration
  67. size_t steps = integrate( harmonic_oscillator ,
  68. x , 0.0 , 10.0 , 0.1 );
  69. //]
  70. //[ integration_class
  71. harm_osc ho(0.15);
  72. steps = integrate( ho ,
  73. x , 0.0 , 10.0 , 0.1 );
  74. //]
  75. //[ integrate_observ
  76. vector<state_type> x_vec;
  77. vector<double> times;
  78. steps = integrate( harmonic_oscillator ,
  79. x , 0.0 , 10.0 , 0.1 ,
  80. push_back_state_and_time( x_vec , times ) );
  81. /* output */
  82. for( size_t i=0; i<=steps; i++ )
  83. {
  84. cout << times[i] << '\t' << x_vec[i][0] << '\t' << x_vec[i][1] << '\n';
  85. }
  86. //]
  87. //[ define_const_stepper
  88. runge_kutta4< state_type > stepper;
  89. integrate_const( stepper , harmonic_oscillator , x , 0.0 , 10.0 , 0.01 );
  90. //]
  91. //[ integrate_const_loop
  92. const double dt = 0.01;
  93. for( double t=0.0 ; t<10.0 ; t+= dt )
  94. stepper.do_step( harmonic_oscillator , x , t , dt );
  95. //]
  96. //[ define_adapt_stepper
  97. typedef runge_kutta_cash_karp54< state_type > error_stepper_type;
  98. //]
  99. //[ integrate_adapt
  100. typedef controlled_runge_kutta< error_stepper_type > controlled_stepper_type;
  101. controlled_stepper_type controlled_stepper;
  102. integrate_adaptive( controlled_stepper , harmonic_oscillator , x , 0.0 , 10.0 , 0.01 );
  103. //]
  104. {
  105. //[integrate_adapt_full
  106. double abs_err = 1.0e-10 , rel_err = 1.0e-6 , a_x = 1.0 , a_dxdt = 1.0;
  107. controlled_stepper_type controlled_stepper(
  108. default_error_checker< double , range_algebra , default_operations >( abs_err , rel_err , a_x , a_dxdt ) );
  109. integrate_adaptive( controlled_stepper , harmonic_oscillator , x , 0.0 , 10.0 , 0.01 );
  110. //]
  111. }
  112. //[integrate_adapt_make_controlled
  113. integrate_adaptive( make_controlled< error_stepper_type >( 1.0e-10 , 1.0e-6 ) ,
  114. harmonic_oscillator , x , 0.0 , 10.0 , 0.01 );
  115. //]
  116. //[integrate_adapt_make_controlled_alternative
  117. integrate_adaptive( make_controlled( 1.0e-10 , 1.0e-6 , error_stepper_type() ) ,
  118. harmonic_oscillator , x , 0.0 , 10.0 , 0.01 );
  119. //]
  120. #ifdef BOOST_NUMERIC_ODEINT_CXX11
  121. //[ define_const_stepper_cpp11
  122. {
  123. runge_kutta4< state_type > stepper;
  124. integrate_const( stepper , []( const state_type &x , state_type &dxdt , double t ) {
  125. dxdt[0] = x[1]; dxdt[1] = -x[0] - gam*x[1]; }
  126. , x , 0.0 , 10.0 , 0.01 );
  127. }
  128. //]
  129. //[ harm_iterator_const_step]
  130. std::for_each( make_const_step_time_iterator_begin( stepper , harmonic_oscillator, x , 0.0 , 0.1 , 10.0 ) ,
  131. make_const_step_time_iterator_end( stepper , harmonic_oscillator, x ) ,
  132. []( std::pair< const state_type & , const double & > x ) {
  133. cout << x.second << " " << x.first[0] << " " << x.first[1] << "\n"; } );
  134. //]
  135. #endif
  136. }