Imported Upstream version 1.57.0
[platform/upstream/boost.git] / libs / numeric / odeint / test / adams_bashforth_moulton.cpp
1 /*
2  [auto_generated]
3  libs/numeric/odeint/test/adams_bashforth_moulton.cpp
4
5  [begin_description]
6  This file tests the use of the Adams-Bashforth-Moulton.
7  [end_description]
8
9  Copyright 2011-2012 Karsten Ahnert
10  Copyright 2011 Mario Mulansky
11
12  Distributed under the Boost Software License, Version 1.0.
13  (See accompanying file LICENSE_1_0.txt or
14  copy at http://www.boost.org/LICENSE_1_0.txt)
15  */
16
17 #define BOOST_TEST_MODULE odeint_adams_bashforth_moulton
18
19 #include <utility>
20
21 #include <boost/test/unit_test.hpp>
22
23 #include <boost/mpl/list.hpp>
24 #include <boost/mpl/size_t.hpp>
25 #include <boost/mpl/range_c.hpp>
26
27 #include <boost/numeric/odeint/stepper/adams_bashforth_moulton.hpp>
28
29 using namespace boost::unit_test;
30 using namespace boost::numeric::odeint;
31
32 typedef double value_type;
33
34 struct lorenz
35 {
36     template< class State , class Deriv , class Value >
37     void operator()( const State &_x , Deriv &_dxdt , const Value &dt ) const
38     {
39         const value_type sigma = 10.0;
40         const value_type R = 28.0;
41         const value_type b = 8.0 / 3.0;
42
43         typename boost::range_iterator< const State >::type x = boost::begin( _x );
44         typename boost::range_iterator< Deriv >::type dxdt = boost::begin( _dxdt );
45
46         dxdt[0] = sigma * ( x[1] - x[0] );
47         dxdt[1] = R * x[0] - x[1] - x[0] * x[2];
48         dxdt[2] = x[0]*x[1] - b * x[2];
49     }
50 };
51
52 BOOST_AUTO_TEST_SUITE( adams_bashforth_moulton_test )
53
54 typedef boost::mpl::range_c< size_t , 1 , 6 > vector_of_steps;
55
56 BOOST_AUTO_TEST_CASE_TEMPLATE( test_init_and_steps , step_type , vector_of_steps )
57 {
58     const static size_t steps = step_type::value;
59     typedef boost::array< value_type , 3 > state_type;
60
61     adams_bashforth_moulton< steps , state_type > stepper;
62     state_type x = {{ 10.0 , 10.0 , 10.0 }};
63     const value_type dt = 0.01;
64     value_type t = 0.0;
65
66     stepper.initialize( lorenz() , x , t , dt );
67     BOOST_CHECK_CLOSE( t , value_type( steps - 1 ) * dt , 1.0e-14 );
68
69     stepper.do_step( lorenz() , x , t , dt );
70 }
71
72
73 BOOST_AUTO_TEST_CASE( test_copying )
74 {
75     typedef boost::array< double , 1 > state_type;
76     typedef adams_bashforth_moulton< 2 , state_type > stepper_type;
77
78     stepper_type s1;
79
80     stepper_type s2( s1 );
81
82     stepper_type s3;
83     s3 = s1;
84  }
85
86
87 BOOST_AUTO_TEST_CASE( test_instantiation )
88 {
89     typedef boost::array< double , 3 > state_type;
90     adams_bashforth_moulton< 1 , state_type > s1;
91     adams_bashforth_moulton< 2 , state_type > s2;
92     adams_bashforth_moulton< 3 , state_type > s3;
93     adams_bashforth_moulton< 4 , state_type > s4;
94     adams_bashforth_moulton< 5 , state_type > s5;
95     adams_bashforth_moulton< 6 , state_type > s6;
96     adams_bashforth_moulton< 7 , state_type > s7;
97     adams_bashforth_moulton< 8 , state_type > s8;
98
99     state_type x = {{ 10.0 , 10.0 , 10.0 }};
100     value_type t = 0.0 , dt = 0.01;
101     s1.do_step( lorenz() , x , t , dt );
102     s2.do_step( lorenz() , x , t , dt );
103     s3.do_step( lorenz() , x , t , dt );
104     s4.do_step( lorenz() , x , t , dt );
105     s5.do_step( lorenz() , x , t , dt );
106     s6.do_step( lorenz() , x , t , dt );
107 //  s7.do_step( lorenz() , x , t , dt );
108 //  s8.do_step( lorenz() , x , t , dt );
109 }
110
111
112 BOOST_AUTO_TEST_SUITE_END()