9
3

lorenz.cpp 4.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147
  1. /* Boost lorenz.cpp test file
  2. Copyright 2012 Karsten Ahnert
  3. Copyright 2012 Mario Mulansky
  4. This file tests the odeint library with the vexcl 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_vexcl
  10. #include <vector>
  11. #include <iostream>
  12. #include <vexcl/vexcl.hpp>
  13. #include <boost/numeric/odeint/stepper/runge_kutta4.hpp>
  14. #include <boost/numeric/odeint/stepper/runge_kutta_cash_karp54.hpp>
  15. #include <boost/numeric/odeint/stepper/runge_kutta_dopri5.hpp>
  16. #include <boost/numeric/odeint/stepper/controlled_runge_kutta.hpp>
  17. #include <boost/numeric/odeint/stepper/dense_output_runge_kutta.hpp>
  18. #include <boost/numeric/odeint/algebra/vector_space_algebra.hpp>
  19. #include <boost/numeric/odeint/integrate/integrate_const.hpp>
  20. #include <boost/numeric/odeint/external/vexcl/vexcl.hpp>
  21. #include <boost/test/unit_test.hpp>
  22. using namespace boost::unit_test;
  23. namespace odeint = boost::numeric::odeint;
  24. typedef vex::vector< double > vector_type;
  25. typedef vex::multivector< double, 3 > state_type;
  26. const double sigma = 10.0;
  27. const double b = 8.0 / 3.0;
  28. struct sys_func
  29. {
  30. const vector_type &R;
  31. sys_func( const vector_type &_R ) : R( _R ) { }
  32. void operator()( const state_type &x , state_type &dxdt , double t ) const
  33. {
  34. dxdt(0) = -sigma * ( x(0) - x(1) );
  35. dxdt(1) = R * x(0) - x(1) - x(0) * x(2);
  36. dxdt(2) = - b * x(2) + x(0) * x(1);
  37. }
  38. };
  39. const size_t n = 1024 ;
  40. const double dt = 0.01;
  41. const double t_max = 1.0;
  42. BOOST_AUTO_TEST_CASE( integrate_const_rk4 )
  43. {
  44. vex::Context ctx( vex::Filter::Type(CL_DEVICE_TYPE_GPU) );
  45. double Rmin = 0.1 , Rmax = 50.0 , dR = ( Rmax - Rmin ) / double( n - 1 );
  46. std::vector<double> x( n * 3 ) , r( n );
  47. for( size_t i=0 ; i<n ; ++i ) r[i] = Rmin + dR * double( i );
  48. state_type X( ctx.queue() , n );
  49. X(0) = 10.0;
  50. X(1) = 10.0;
  51. X(2) = 10.0;
  52. vector_type R( ctx.queue() , r );
  53. odeint::runge_kutta4<
  54. state_type , double , state_type , double ,
  55. odeint::vector_space_algebra , odeint::default_operations
  56. > stepper;
  57. odeint::integrate_const( stepper , sys_func( R ) , X , 0.0 , t_max , dt );
  58. std::vector< double > res( 3 * n );
  59. vex::copy( X(0) , res );
  60. std::cout << res[0] << std::endl;
  61. }
  62. BOOST_AUTO_TEST_CASE( integrate_const_controlled_rk54 )
  63. {
  64. vex::Context ctx( vex::Filter::Type(CL_DEVICE_TYPE_GPU) );
  65. double Rmin = 0.1 , Rmax = 50.0 , dR = ( Rmax - Rmin ) / double( n - 1 );
  66. std::vector<double> x( n * 3 ) , r( n );
  67. for( size_t i=0 ; i<n ; ++i ) r[i] = Rmin + dR * double( i );
  68. state_type X( ctx.queue() , n );
  69. X(0) = 10.0;
  70. X(1) = 10.0;
  71. X(2) = 10.0;
  72. vector_type R( ctx.queue() , r );
  73. typedef odeint::runge_kutta_cash_karp54<
  74. state_type , double , state_type , double ,
  75. odeint::vector_space_algebra , odeint::default_operations
  76. > stepper_type;
  77. typedef odeint::controlled_runge_kutta< stepper_type > controlled_stepper_type;
  78. odeint::integrate_const( controlled_stepper_type() , sys_func( R ) , X , 0.0 , t_max , dt );
  79. std::vector< double > res( 3 * n );
  80. vex::copy( X(0) , res );
  81. std::cout << res[0] << std::endl;
  82. }
  83. BOOST_AUTO_TEST_CASE( integrate_const_dense_output_dopri5 )
  84. {
  85. vex::Context ctx( vex::Filter::Type(CL_DEVICE_TYPE_GPU) );
  86. double Rmin = 0.1 , Rmax = 50.0 , dR = ( Rmax - Rmin ) / double( n - 1 );
  87. std::vector<double> x( n * 3 ) , r( n );
  88. for( size_t i=0 ; i<n ; ++i ) r[i] = Rmin + dR * double( i );
  89. state_type X( ctx.queue() , n );
  90. X(0) = 10.0;
  91. X(1) = 10.0;
  92. X(2) = 10.0;
  93. vector_type R( ctx.queue() , r );
  94. typedef odeint::runge_kutta_dopri5<
  95. state_type , double , state_type , double ,
  96. odeint::vector_space_algebra , odeint::default_operations
  97. > stepper_type;
  98. typedef odeint::controlled_runge_kutta< stepper_type > controlled_stepper_type;
  99. typedef odeint::dense_output_runge_kutta< controlled_stepper_type > dense_output_stepper_type;
  100. odeint::integrate_const( dense_output_stepper_type() , sys_func( R ) , X , 0.0 , t_max , dt );
  101. std::vector< double > res( 3 * n );
  102. vex::copy( X(0) , res );
  103. std::cout << res[0] << std::endl;
  104. }