• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
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 
14 #include <iostream>
15 #include <cmath>
16 
17 #include <boost/array.hpp>
18 #include <boost/numeric/odeint.hpp>
19 
20 using namespace boost::numeric::odeint;
21 
22 const int Steps = 4;
23 
24 typedef double value_type;
25 
26 typedef boost::array< double , 2 > state_type;
27 
28 typedef runge_kutta_fehlberg78<state_type> initializing_stepper_type;
29 typedef adams_bashforth_moulton< Steps , state_type > stepper_type;
30 //typedef adams_bashforth< Steps , state_type > stepper_type;
31 
32 // harmonic oscillator, analytic solution x[0] = sin( t )
33 struct osc
34 {
operator ()osc35     void operator()( const state_type &x , state_type &dxdt , const double t ) const
36     {
37         dxdt[0] = x[1];
38         dxdt[1] = -x[0];
39     }
40 };
41 
main()42 int main()
43 {
44     stepper_type stepper;
45     initializing_stepper_type init_stepper;
46     const int o = stepper.order()+1; //order of the error is order of approximation + 1
47 
48     const state_type x0 = {{ 0.0 , 1.0 }};
49     state_type x1 = x0;
50     double t = 0.0;
51     double dt = 0.25;
52     // initialization, does a number of steps already to fill internal buffer, t is increased
53     // we use the rk78 as initializing stepper
54     stepper.initialize( boost::ref(init_stepper) , osc() , x1 , t , dt );
55     // do a number of steps to fill the buffer with results from adams bashforth
56     for( size_t n=0 ; n < stepper.steps ; ++n )
57     {
58         stepper.do_step( osc() , x1 , t , dt );
59         t += dt;
60     }
61     double A = std::sqrt( x1[0]*x1[0] + x1[1]*x1[1] );
62     double phi = std::asin(x1[0]/A) - t;
63     // now we do the actual step
64     stepper.do_step( osc() , x1 , t , dt );
65     // only examine the error of the adams-bashforth-moulton step, not the initialization
66     const double f = 2.0 * std::abs( A*sin(t+dt+phi) - x1[0] ) / std::pow( dt , o ); // upper bound
67 
68     std::cout << "# " << o << " , " << f << std::endl;
69 
70     /* as long as we have errors above machine precision */
71     while( f*std::pow( dt , o ) > 1E-16 )
72     {
73         x1 = x0;
74         t = 0.0;
75         stepper.initialize( boost::ref(init_stepper) , osc() , x1 , t , dt );
76         A = std::sqrt( x1[0]*x1[0] + x1[1]*x1[1] );
77         phi = std::asin(x1[0]/A) - t;
78         // now we do the actual step
79         stepper.do_step( osc() , x1 , t , dt );
80         // only examine the error of the adams-bashforth-moulton step, not the initialization
81         std::cout << dt << '\t' <<  std::abs( A*sin(t+dt+phi) - x1[0] ) << std::endl;
82         dt *= 0.5;
83     }
84 }
85