• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
1 /*
2 
3  [begin_description]
4  Test case for issue 189:
5  Controlled Rosenbrock stepper fails to increase step size
6  [end_description]
7 
8  Copyright 2016 Karsten Ahnert
9  Copyright 2016 Mario Mulansky
10 
11  Distributed under the Boost Software License, Version 1.0.
12  (See accompanying file LICENSE_1_0.txt or
13  copy at http://www.boost.org/LICENSE_1_0.txt)
14  */
15 
16 #define BOOST_TEST_MODULE odeint_regression_189
17 
18 #include <boost/numeric/odeint.hpp>
19 
20 #include <boost/test/unit_test.hpp>
21 #include <boost/phoenix/core.hpp>
22 #include <boost/phoenix/operator.hpp>
23 
24 using namespace boost::numeric::odeint;
25 namespace phoenix = boost::phoenix;
26 
27 typedef boost::numeric::ublas::vector< double > vector_type;
28 typedef boost::numeric::ublas::matrix< double > matrix_type;
29 
30 struct stiff_system
31 {
operator ()stiff_system32     void operator()( const vector_type &x , vector_type &dxdt , double /* t */ )
33     {
34         dxdt[ 0 ] = -101.0 * x[ 0 ] - 100.0 * x[ 1 ];
35         dxdt[ 1 ] = x[ 0 ];
36     }
37 };
38 
39 struct stiff_system_jacobi
40 {
operator ()stiff_system_jacobi41     void operator()( const vector_type & /* x */ , matrix_type &J , const double & /* t */ , vector_type &dfdt )
42     {
43         J( 0 , 0 ) = -101.0;
44         J( 0 , 1 ) = -100.0;
45         J( 1 , 0 ) = 1.0;
46         J( 1 , 1 ) = 0.0;
47         dfdt[0] = 0.0;
48         dfdt[1] = 0.0;
49     }
50 };
51 
52 
BOOST_AUTO_TEST_CASE(regression_189)53 BOOST_AUTO_TEST_CASE( regression_189 )
54 {
55     vector_type x( 2 , 1.0 );
56 
57     size_t num_of_steps = integrate_const( make_dense_output< rosenbrock4< double > >( 1.0e-6 , 1.0e-6 ) ,
58             std::make_pair( stiff_system() , stiff_system_jacobi() ) ,
59             x , 0.0 , 50.0 , 0.01 ,
60             std::cout << phoenix::arg_names::arg2 << " " << phoenix::arg_names::arg1[0] << "\n" );
61     // regression: number of steps should be 74
62     size_t num_of_steps_expected = 74;
63     BOOST_CHECK_EQUAL( num_of_steps , num_of_steps_expected );
64 
65     vector_type x2( 2 , 1.0 );
66 
67     size_t num_of_steps2 = integrate_const( make_dense_output< runge_kutta_dopri5< vector_type > >( 1.0e-6 , 1.0e-6 ) ,
68             stiff_system() , x2 , 0.0 , 50.0 , 0.01 ,
69             std::cout << phoenix::arg_names::arg2 << " " << phoenix::arg_names::arg1[0] << "\n" );
70     num_of_steps_expected = 1531;
71     BOOST_CHECK_EQUAL( num_of_steps2 , num_of_steps_expected );
72 }
73