lorenz_point.cpp 2.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119
  1. /*
  2. * Copyright 2011-2013 Mario Mulansky
  3. * Copyright 2012 Karsten Ahnert
  4. *
  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. * Example for the lorenz system with a 3D point type
  10. */
  11. #include <iostream>
  12. #include <cmath>
  13. #include <boost/operators.hpp>
  14. #include <boost/numeric/odeint.hpp>
  15. //[point3D
  16. class point3D :
  17. boost::additive1< point3D ,
  18. boost::additive2< point3D , double ,
  19. boost::multiplicative2< point3D , double > > >
  20. {
  21. public:
  22. double x , y , z;
  23. point3D()
  24. : x( 0.0 ) , y( 0.0 ) , z( 0.0 )
  25. { }
  26. point3D( const double val )
  27. : x( val ) , y( val ) , z( val )
  28. { }
  29. point3D( const double _x , const double _y , const double _z )
  30. : x( _x ) , y( _y ) , z( _z )
  31. { }
  32. point3D& operator+=( const point3D &p )
  33. {
  34. x += p.x; y += p.y; z += p.z;
  35. return *this;
  36. }
  37. point3D& operator*=( const double a )
  38. {
  39. x *= a; y *= a; z *= a;
  40. return *this;
  41. }
  42. };
  43. //]
  44. //[point3D_abs_div
  45. // only required for steppers with error control
  46. point3D operator/( const point3D &p1 , const point3D &p2 )
  47. {
  48. return point3D( p1.x/p2.x , p1.y/p2.y , p1.z/p2.z );
  49. }
  50. point3D abs( const point3D &p )
  51. {
  52. return point3D( std::abs(p.x) , std::abs(p.y) , std::abs(p.z) );
  53. }
  54. //]
  55. //[point3D_norm
  56. // also only for steppers with error control
  57. namespace boost { namespace numeric { namespace odeint {
  58. template<>
  59. struct vector_space_norm_inf< point3D >
  60. {
  61. typedef double result_type;
  62. double operator()( const point3D &p ) const
  63. {
  64. using std::max;
  65. using std::abs;
  66. return max( max( abs( p.x ) , abs( p.y ) ) , abs( p.z ) );
  67. }
  68. };
  69. } } }
  70. //]
  71. std::ostream& operator<<( std::ostream &out , const point3D &p )
  72. {
  73. out << p.x << " " << p.y << " " << p.z;
  74. return out;
  75. }
  76. //[point3D_main
  77. const double sigma = 10.0;
  78. const double R = 28.0;
  79. const double b = 8.0 / 3.0;
  80. void lorenz( const point3D &x , point3D &dxdt , const double t )
  81. {
  82. dxdt.x = sigma * ( x.y - x.x );
  83. dxdt.y = R * x.x - x.y - x.x * x.z;
  84. dxdt.z = -b * x.z + x.x * x.y;
  85. }
  86. using namespace boost::numeric::odeint;
  87. int main()
  88. {
  89. point3D x( 10.0 , 5.0 , 5.0 );
  90. // point type defines it's own operators -> use vector_space_algebra !
  91. typedef runge_kutta_dopri5< point3D , double , point3D ,
  92. double , vector_space_algebra > stepper;
  93. int steps = integrate_adaptive( make_controlled<stepper>( 1E-10 , 1E-10 ) , lorenz , x ,
  94. 0.0 , 10.0 , 0.1 );
  95. std::cout << x << std::endl;
  96. std::cout << "steps: " << steps << std::endl;
  97. }
  98. //]