abm_precision.cpp 2.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384
  1. /*
  2. * abm_precision.cpp
  3. *
  4. * example to check the order of the multi-step methods
  5. *
  6. * Copyright 2009-2013 Karsten Ahnert
  7. * Copyright 2009-2013 Mario Mulansky
  8. *
  9. * Distributed under the Boost Software License, Version 1.0.
  10. * (See accompanying file LICENSE_1_0.txt or
  11. * copy at http://www.boost.org/LICENSE_1_0.txt)
  12. */
  13. #include <iostream>
  14. #include <cmath>
  15. #include <boost/array.hpp>
  16. #include <boost/numeric/odeint.hpp>
  17. using namespace boost::numeric::odeint;
  18. const int Steps = 4;
  19. typedef double value_type;
  20. typedef boost::array< double , 2 > state_type;
  21. typedef runge_kutta_fehlberg78<state_type> initializing_stepper_type;
  22. typedef adams_bashforth_moulton< Steps , state_type > stepper_type;
  23. //typedef adams_bashforth< Steps , state_type > stepper_type;
  24. // harmonic oscillator, analytic solution x[0] = sin( t )
  25. struct osc
  26. {
  27. void operator()( const state_type &x , state_type &dxdt , const double t ) const
  28. {
  29. dxdt[0] = x[1];
  30. dxdt[1] = -x[0];
  31. }
  32. };
  33. int main()
  34. {
  35. stepper_type stepper;
  36. initializing_stepper_type init_stepper;
  37. const int o = stepper.order()+1; //order of the error is order of approximation + 1
  38. const state_type x0 = {{ 0.0 , 1.0 }};
  39. state_type x1 = x0;
  40. double t = 0.0;
  41. double dt = 0.25;
  42. // initialization, does a number of steps already to fill internal buffer, t is increased
  43. // we use the rk78 as initializing stepper
  44. stepper.initialize( boost::ref(init_stepper) , osc() , x1 , t , dt );
  45. // do a number of steps to fill the buffer with results from adams bashforth
  46. for( size_t n=0 ; n < stepper.steps ; ++n )
  47. {
  48. stepper.do_step( osc() , x1 , t , dt );
  49. t += dt;
  50. }
  51. double A = std::sqrt( x1[0]*x1[0] + x1[1]*x1[1] );
  52. double phi = std::asin(x1[0]/A) - t;
  53. // now we do the actual step
  54. stepper.do_step( osc() , x1 , t , dt );
  55. // only examine the error of the adams-bashforth-moulton step, not the initialization
  56. const double f = 2.0 * std::abs( A*sin(t+dt+phi) - x1[0] ) / std::pow( dt , o ); // upper bound
  57. std::cout << "# " << o << " , " << f << std::endl;
  58. /* as long as we have errors above machine precision */
  59. while( f*std::pow( dt , o ) > 1E-16 )
  60. {
  61. x1 = x0;
  62. t = 0.0;
  63. stepper.initialize( boost::ref(init_stepper) , osc() , x1 , t , dt );
  64. A = std::sqrt( x1[0]*x1[0] + x1[1]*x1[1] );
  65. phi = std::asin(x1[0]/A) - t;
  66. // now we do the actual step
  67. stepper.do_step( osc() , x1 , t , dt );
  68. // only examine the error of the adams-bashforth-moulton step, not the initialization
  69. std::cout << dt << '\t' << std::abs( A*sin(t+dt+phi) - x1[0] ) << std::endl;
  70. dt *= 0.5;
  71. }
  72. }