gmp_integrate.cpp 5.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169
  1. /* Boost check_gmp.cpp test file
  2. Copyright 2010-2012 Mario Mulansky
  3. Copyright 2011-2012 Karsten Ahnert
  4. This file tests the odeint library with the gmp arbitrary precision types
  5. Distributed under the Boost Software License, Version 1.0.
  6. (See accompanying file LICENSE_1_0.txt or
  7. copy at http://www.boost.org/LICENSE_1_0.txt)
  8. */
  9. #define BOOST_TEST_MODULE odeint_gmp
  10. #include <gmpxx.h>
  11. #include <boost/test/unit_test.hpp>
  12. #include <boost/array.hpp>
  13. #include <boost/mpl/vector.hpp>
  14. #include <boost/numeric/odeint.hpp>
  15. //#include <boost/numeric/odeint/algebra/vector_space_algebra.hpp>
  16. using namespace boost::unit_test;
  17. using namespace boost::numeric::odeint;
  18. namespace mpl = boost::mpl;
  19. const int precision = 1024;
  20. typedef mpf_class value_type;
  21. typedef mpf_class state_type;
  22. //provide min, max and pow functions for mpf types - required for controlled steppers
  23. value_type min( const value_type a , const value_type b )
  24. {
  25. if( a<b ) return a;
  26. else return b;
  27. }
  28. value_type max( const value_type a , const value_type b )
  29. {
  30. if( a>b ) return a;
  31. else return b;
  32. }
  33. value_type pow( const value_type a , const value_type b )
  34. {
  35. // do the calculation in double precision...
  36. return value_type( std::pow( a.get_d() , b.get_d() ) );
  37. }
  38. //provide vector_space reduce:
  39. namespace boost { namespace numeric { namespace odeint {
  40. template<>
  41. struct vector_space_reduce< state_type >
  42. {
  43. template< class Op >
  44. state_type operator()( state_type x , Op op , state_type init ) const
  45. {
  46. init = op( init , x );
  47. return init;
  48. }
  49. };
  50. } } }
  51. void constant_system( const state_type &x , state_type &dxdt , value_type t )
  52. {
  53. dxdt = value_type( 1.0 , precision );
  54. }
  55. /* check runge kutta stepers */
  56. typedef mpl::vector<
  57. euler< state_type , value_type , state_type , value_type , vector_space_algebra > ,
  58. modified_midpoint< state_type , value_type , state_type , value_type , vector_space_algebra > ,
  59. runge_kutta4< state_type , value_type , state_type , value_type , vector_space_algebra > ,
  60. runge_kutta4_classic< state_type , value_type , state_type , value_type , vector_space_algebra > ,
  61. runge_kutta_cash_karp54_classic< state_type , value_type , state_type , value_type , vector_space_algebra > ,
  62. runge_kutta_cash_karp54< state_type , value_type , state_type , value_type , vector_space_algebra > ,
  63. runge_kutta_dopri5< state_type , value_type , state_type , value_type , vector_space_algebra > ,
  64. runge_kutta_fehlberg78< state_type , value_type , state_type , value_type , vector_space_algebra >
  65. > stepper_types;
  66. template< class Stepper >
  67. struct perform_integrate_const_test {
  68. void operator()( void )
  69. {
  70. /* We have to specify the desired precision in advance! */
  71. mpf_set_default_prec( precision );
  72. mpf_t eps_ , unity;
  73. mpf_init( eps_ ); mpf_init( unity );
  74. mpf_set_d( unity , 1.0 );
  75. mpf_div_2exp( eps_ , unity , precision-1 ); // 2^(-precision+1) : smallest number that can be represented with used precision
  76. value_type eps( eps_ );
  77. Stepper stepper;
  78. state_type x;
  79. x = 0.0;
  80. value_type t0( 0.0 );
  81. value_type tend( 1.0 );
  82. value_type dt(0.1);
  83. integrate_const( stepper , constant_system , x , t0 , tend , dt );
  84. x = 0.0;
  85. t0 = 0.0;
  86. dt = 0.1;
  87. size_t steps = 10;
  88. integrate_n_steps( stepper , constant_system , x , t0 , dt , steps );
  89. BOOST_CHECK_MESSAGE( abs( x - 10*dt ) < eps , x );
  90. }
  91. };
  92. BOOST_AUTO_TEST_CASE_TEMPLATE( integrate_const_test , Stepper , stepper_types )
  93. {
  94. perform_integrate_const_test< Stepper > tester;
  95. tester();
  96. }
  97. typedef mpl::vector<
  98. controlled_runge_kutta< runge_kutta_cash_karp54_classic< state_type , value_type , state_type , value_type , vector_space_algebra > > ,
  99. controlled_runge_kutta< runge_kutta_dopri5< state_type , value_type , state_type , value_type , vector_space_algebra > > ,
  100. controlled_runge_kutta< runge_kutta_fehlberg78< state_type , value_type , state_type , value_type , vector_space_algebra > > ,
  101. bulirsch_stoer< state_type , value_type , state_type , value_type , vector_space_algebra >
  102. > controlled_stepper_types;
  103. template< class Stepper >
  104. struct perform_integrate_adaptive_test {
  105. void operator()( void )
  106. {
  107. mpf_set_default_prec( precision );
  108. mpf_t eps_ , unity;
  109. mpf_init( eps_ ); mpf_init( unity );
  110. mpf_set_d( unity , 1.0 );
  111. mpf_div_2exp( eps_ , unity , precision-1 ); // 2^(-precision+1) : smallest number that can be represented with used precision
  112. value_type eps( eps_ );
  113. Stepper stepper;
  114. state_type x;
  115. x = 0.0;
  116. value_type t0( 0.0 );
  117. value_type tend( 1.0 );
  118. value_type dt(0.1);
  119. integrate_adaptive( stepper , constant_system , x , t0 , tend , dt );
  120. BOOST_CHECK_MESSAGE( abs( x - tend ) < eps , x - 0.1 );
  121. }
  122. };
  123. BOOST_AUTO_TEST_CASE_TEMPLATE( integrate_adaptive__test , Stepper , controlled_stepper_types )
  124. {
  125. perform_integrate_adaptive_test< Stepper > tester;
  126. tester();
  127. }