Brian Silverman | 7c33ab2 | 2018-08-04 17:14:51 -0700 | [diff] [blame^] | 1 | /* |
| 2 | [auto_generated] |
| 3 | libs/numeric/odeint/test/integrate.cpp |
| 4 | |
| 5 | [begin_description] |
| 6 | This file tests the integrate function and its variants. |
| 7 | [end_description] |
| 8 | |
| 9 | Copyright 2011-2012 Mario Mulansky |
| 10 | Copyright 2011-2012 Karsten Ahnert |
| 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 | |
| 18 | #define BOOST_TEST_MODULE odeint_integrate_functions |
| 19 | |
| 20 | #include <vector> |
| 21 | #include <cmath> |
| 22 | #include <iostream> |
| 23 | |
| 24 | #include <boost/numeric/odeint/config.hpp> |
| 25 | |
| 26 | #include <boost/array.hpp> |
| 27 | #include <boost/ref.hpp> |
| 28 | #include <boost/iterator/counting_iterator.hpp> |
| 29 | |
| 30 | #include <boost/test/unit_test.hpp> |
| 31 | |
| 32 | #include <boost/mpl/vector.hpp> |
| 33 | |
| 34 | // nearly everything from odeint is used in these tests |
| 35 | #ifndef ODEINT_INTEGRATE_ITERATOR |
| 36 | #include <boost/numeric/odeint/integrate/integrate_const.hpp> |
| 37 | #include <boost/numeric/odeint/integrate/integrate_adaptive.hpp> |
| 38 | #include <boost/numeric/odeint/integrate/integrate_times.hpp> |
| 39 | #include <boost/numeric/odeint/integrate/integrate_n_steps.hpp> |
| 40 | #else |
| 41 | #include <boost/numeric/odeint/iterator/integrate/integrate_const.hpp> |
| 42 | #include <boost/numeric/odeint/iterator/integrate/integrate_adaptive.hpp> |
| 43 | #include <boost/numeric/odeint/iterator/integrate/integrate_times.hpp> |
| 44 | #include <boost/numeric/odeint/iterator/integrate/integrate_n_steps.hpp> |
| 45 | #endif |
| 46 | #include <boost/numeric/odeint/stepper/euler.hpp> |
| 47 | #include <boost/numeric/odeint/stepper/modified_midpoint.hpp> |
| 48 | #include <boost/numeric/odeint/stepper/runge_kutta4.hpp> |
| 49 | #include <boost/numeric/odeint/stepper/runge_kutta_cash_karp54.hpp> |
| 50 | #include <boost/numeric/odeint/stepper/runge_kutta_dopri5.hpp> |
| 51 | #include <boost/numeric/odeint/stepper/runge_kutta_fehlberg78.hpp> |
| 52 | #include <boost/numeric/odeint/stepper/controlled_runge_kutta.hpp> |
| 53 | #include <boost/numeric/odeint/stepper/bulirsch_stoer.hpp> |
| 54 | #include <boost/numeric/odeint/stepper/dense_output_runge_kutta.hpp> |
| 55 | #include <boost/numeric/odeint/stepper/bulirsch_stoer_dense_out.hpp> |
| 56 | |
| 57 | #include <boost/numeric/odeint/util/detail/less_with_sign.hpp> |
| 58 | |
| 59 | using namespace boost::unit_test; |
| 60 | using namespace boost::numeric::odeint; |
| 61 | namespace mpl = boost::mpl; |
| 62 | |
| 63 | |
| 64 | typedef double value_type; |
| 65 | typedef std::vector< value_type > state_type; |
| 66 | |
| 67 | void lorenz( const state_type &x , state_type &dxdt , const value_type t ) |
| 68 | { |
| 69 | //const value_type sigma( 10.0 ); |
| 70 | const value_type R( 28.0 ); |
| 71 | const value_type b( value_type( 8.0 ) / value_type( 3.0 ) ); |
| 72 | |
| 73 | // first component trivial |
| 74 | dxdt[0] = 1.0; //sigma * ( x[1] - x[0] ); |
| 75 | dxdt[1] = R * x[0] - x[1] - x[0] * x[2]; |
| 76 | dxdt[2] = -b * x[2] + x[0] * x[1]; |
| 77 | } |
| 78 | |
| 79 | struct push_back_time |
| 80 | { |
| 81 | std::vector< double >& m_times; |
| 82 | |
| 83 | state_type& m_x; |
| 84 | |
| 85 | push_back_time( std::vector< double > × , state_type &x ) |
| 86 | : m_times( times ) , m_x( x ) { } |
| 87 | |
| 88 | void operator()( const state_type &x , double t ) |
| 89 | { |
| 90 | m_times.push_back( t ); |
| 91 | boost::numeric::odeint::copy( x , m_x ); |
| 92 | } |
| 93 | }; |
| 94 | |
| 95 | template< class Stepper > |
| 96 | struct perform_integrate_const_test |
| 97 | { |
| 98 | void operator()( const value_type t_end , const value_type dt ) |
| 99 | { |
| 100 | // std::cout << "Testing integrate_const with " << typeid( Stepper ).name() << std::endl; |
| 101 | |
| 102 | state_type x( 3 , 10.0 ) , x_end( 3 ); |
| 103 | |
| 104 | std::vector< value_type > times; |
| 105 | |
| 106 | size_t steps = integrate_const( Stepper() , lorenz , x , 0.0 , t_end , |
| 107 | dt , push_back_time( times , x_end ) ); |
| 108 | |
| 109 | std::cout.precision(16); |
| 110 | std::cout << t_end << " (" << dt << "), " << steps << " , " << times.size() << " , " << 10.0+dt*steps << "=" << x_end[0] << std::endl; |
| 111 | |
| 112 | std::cout << static_cast<int>(floor(t_end/dt)) << " , " << t_end/dt << std::endl; |
| 113 | |
| 114 | BOOST_CHECK_EQUAL( static_cast<int>(times.size()) , static_cast<int>(floor(t_end/dt))+1 ); |
| 115 | |
| 116 | for( size_t i=0 ; i<times.size() ; ++i ) |
| 117 | { |
| 118 | //std::cout << i << " , " << times[i] << " , " << static_cast< value_type >(i)*dt << std::endl; |
| 119 | // check if observer was called at times 0,1,2,... |
| 120 | BOOST_CHECK_SMALL( times[i] - static_cast< value_type >(i)*dt , (i+1) * 2E-16 ); |
| 121 | } |
| 122 | |
| 123 | // check first, trivial, component |
| 124 | BOOST_CHECK_SMALL( (10.0 + times[times.size()-1]) - x_end[0] , 1E-6 ); // precision of steppers: 1E-6 |
| 125 | //BOOST_CHECK_EQUAL( x[1] , x_end[1] ); |
| 126 | //BOOST_CHECK_EQUAL( x[2] , x_end[2] ); |
| 127 | } |
| 128 | }; |
| 129 | |
| 130 | template< class Stepper > |
| 131 | struct perform_integrate_adaptive_test |
| 132 | { |
| 133 | void operator()( const value_type t_end = 10.0 , const value_type dt = 0.03 ) |
| 134 | { |
| 135 | // std::cout << "Testing integrate_adaptive with " << typeid( Stepper ).name() << std::endl; |
| 136 | |
| 137 | state_type x( 3 , 10.0 ) , x_end( 3 ); |
| 138 | |
| 139 | std::vector< value_type > times; |
| 140 | |
| 141 | size_t steps = integrate_adaptive( Stepper() , lorenz , x , 0.0 , t_end , |
| 142 | dt , push_back_time( times , x_end ) ); |
| 143 | |
| 144 | // std::cout << t_end << " , " << steps << " , " << times.size() << " , " << dt << " , " << 10.0+t_end << "=" << x_end[0] << std::endl; |
| 145 | |
| 146 | BOOST_CHECK_EQUAL( times.size() , steps+1 ); |
| 147 | |
| 148 | BOOST_CHECK_SMALL( times[0] - 0.0 , 2E-16 ); |
| 149 | BOOST_CHECK_SMALL( times[times.size()-1] - t_end , times.size() * 2E-16 ); |
| 150 | |
| 151 | // check first, trivial, component |
| 152 | BOOST_CHECK_SMALL( (10.0 + t_end) - x_end[0] , 1E-6 ); // precision of steppers: 1E-6 |
| 153 | // BOOST_CHECK_EQUAL( x[1] , x_end[1] ); |
| 154 | // BOOST_CHECK_EQUAL( x[2] , x_end[2] ); |
| 155 | } |
| 156 | }; |
| 157 | |
| 158 | |
| 159 | template< class Stepper > |
| 160 | struct perform_integrate_times_test |
| 161 | { |
| 162 | void operator()( const int n = 10 , const int dn=1 , const value_type dt = 0.03 ) |
| 163 | { |
| 164 | // std::cout << "Testing integrate_times with " << typeid( Stepper ).name() << std::endl; |
| 165 | |
| 166 | state_type x( 3 ) , x_end( 3 ); |
| 167 | x[0] = x[1] = x[2] = 10.0; |
| 168 | |
| 169 | std::vector< double > times; |
| 170 | |
| 171 | std::vector< double > obs_times( abs(n) ); |
| 172 | for( int i=0 ; boost::numeric::odeint::detail::less_with_sign( static_cast<double>(i) , |
| 173 | static_cast<double>(obs_times.size()) , |
| 174 | dt ) ; i+=dn ) |
| 175 | { |
| 176 | obs_times[i] = i; |
| 177 | } |
| 178 | // simple stepper |
| 179 | integrate_times( Stepper() , lorenz , x , obs_times.begin() , obs_times.end() , |
| 180 | dt , push_back_time( times , x_end ) ); |
| 181 | |
| 182 | BOOST_CHECK_EQUAL( static_cast<int>(times.size()) , abs(n) ); |
| 183 | |
| 184 | for( size_t i=0 ; i<times.size() ; ++i ) |
| 185 | // check if observer was called at times 0,1,2,... |
| 186 | BOOST_CHECK_EQUAL( times[i] , static_cast<double>(i) ); |
| 187 | |
| 188 | // check first, trivial, component |
| 189 | BOOST_CHECK_SMALL( (10.0 + 1.0*times[times.size()-1]) - x_end[0] , 1E-6 ); // precision of steppers: 1E-6 |
| 190 | // BOOST_CHECK_EQUAL( x[1] , x_end[1] ); |
| 191 | // BOOST_CHECK_EQUAL( x[2] , x_end[2] ); |
| 192 | } |
| 193 | }; |
| 194 | |
| 195 | template< class Stepper > |
| 196 | struct perform_integrate_n_steps_test |
| 197 | { |
| 198 | void operator()( const int n = 200 , const value_type dt = 0.01 ) |
| 199 | { |
| 200 | // std::cout << "Testing integrate_n_steps with " << typeid( Stepper ).name() << ". "; |
| 201 | // std::cout << "dt=" << dt << std::endl; |
| 202 | |
| 203 | state_type x( 3 ) , x_end( 3 ); |
| 204 | x[0] = x[1] = x[2] = 10.0; |
| 205 | |
| 206 | std::vector< double > times; |
| 207 | |
| 208 | // simple stepper |
| 209 | value_type end_time = integrate_n_steps( Stepper() , lorenz , x , 0.0 , dt , n , push_back_time( times , x_end ) ); |
| 210 | |
| 211 | BOOST_CHECK_SMALL( end_time - n*dt , 2E-16 ); |
| 212 | BOOST_CHECK_EQUAL( static_cast<int>(times.size()) , n+1 ); |
| 213 | |
| 214 | for( size_t i=0 ; i<times.size() ; ++i ) |
| 215 | { |
| 216 | // check if observer was called at times 0,1,2,... |
| 217 | BOOST_CHECK_SMALL(times[i] - static_cast< value_type >(i) * dt, 2E-16); |
| 218 | } |
| 219 | // check first, trivial, component |
| 220 | BOOST_CHECK_SMALL( (10.0 + end_time) - x_end[0] , 1E-6 ); // precision of steppers: 1E-6 |
| 221 | // BOOST_CHECK_EQUAL( x[1] , x_end[1] ); |
| 222 | // BOOST_CHECK_EQUAL( x[2] , x_end[2] ); |
| 223 | |
| 224 | } |
| 225 | }; |
| 226 | |
| 227 | |
| 228 | |
| 229 | class stepper_methods : public mpl::vector< |
| 230 | euler< state_type > , |
| 231 | modified_midpoint< state_type > , |
| 232 | runge_kutta4< state_type > , |
| 233 | runge_kutta_cash_karp54< state_type > , |
| 234 | runge_kutta_dopri5< state_type > , |
| 235 | runge_kutta_fehlberg78< state_type > , |
| 236 | controlled_runge_kutta< runge_kutta_cash_karp54< state_type > > , |
| 237 | controlled_runge_kutta< runge_kutta_dopri5< state_type > > , |
| 238 | controlled_runge_kutta< runge_kutta_fehlberg78< state_type > > , |
| 239 | bulirsch_stoer< state_type > , |
| 240 | dense_output_runge_kutta< controlled_runge_kutta< runge_kutta_dopri5< state_type > > > |
| 241 | //bulirsch_stoer_dense_out< state_type > |
| 242 | > { }; |
| 243 | |
| 244 | |
| 245 | |
| 246 | BOOST_AUTO_TEST_SUITE( integrate_test ) |
| 247 | |
| 248 | BOOST_AUTO_TEST_CASE_TEMPLATE( integrate_const_test_case , Stepper, stepper_methods ) |
| 249 | { |
| 250 | perform_integrate_const_test< Stepper > tester; |
| 251 | tester( 1.005 , 0.01 ); |
| 252 | tester( 1.0 , 0.01 ); |
| 253 | tester( 1.1 , 0.01 ); |
| 254 | tester( -1.005 , -0.01 ); |
| 255 | } |
| 256 | |
| 257 | |
| 258 | BOOST_AUTO_TEST_CASE_TEMPLATE( integrate_adaptive_test_case , Stepper, stepper_methods ) |
| 259 | { |
| 260 | perform_integrate_adaptive_test< Stepper > tester; |
| 261 | tester( 1.005 , 0.01 ); |
| 262 | tester( 1.0 , 0.01 ); |
| 263 | tester( 1.1 , 0.01 ); |
| 264 | tester( -1.005 , -0.01 ); |
| 265 | } |
| 266 | |
| 267 | |
| 268 | BOOST_AUTO_TEST_CASE_TEMPLATE( integrate_times_test_case , Stepper, stepper_methods ) |
| 269 | { |
| 270 | perform_integrate_times_test< Stepper > tester; |
| 271 | tester(); |
| 272 | //tester( -10 , -0.01 ); |
| 273 | } |
| 274 | |
| 275 | BOOST_AUTO_TEST_CASE_TEMPLATE( integrate_n_steps_test_case , Stepper, stepper_methods ) |
| 276 | { |
| 277 | perform_integrate_n_steps_test< Stepper > tester; |
| 278 | tester(); |
| 279 | tester( 200 , 0.01 ); |
| 280 | tester( 200 , 0.01 ); |
| 281 | tester( 200 , 0.01 ); |
| 282 | tester( 200 , -0.01 ); |
| 283 | } |
| 284 | |
| 285 | BOOST_AUTO_TEST_SUITE_END() |