lorenz_parameters.cu 10 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296
  1. /*
  2. Copyright 2011-2012 Karsten Ahnert
  3. Copyright 2011-2013 Mario Mulansky
  4. Distributed under the Boost Software License, Version 1.0.
  5. (See accompanying file LICENSE_1_0.txt or
  6. copy at http://www.boost.org/LICENSE_1_0.txt)
  7. */
  8. #include <iostream>
  9. #include <cmath>
  10. #include <utility>
  11. #include <thrust/device_vector.h>
  12. #include <thrust/reduce.h>
  13. #include <thrust/functional.h>
  14. #include <boost/numeric/odeint.hpp>
  15. #include <boost/numeric/odeint/external/thrust/thrust.hpp>
  16. #include <boost/random/mersenne_twister.hpp>
  17. #include <boost/random/uniform_real.hpp>
  18. #include <boost/random/variate_generator.hpp>
  19. using namespace std;
  20. using namespace boost::numeric::odeint;
  21. //change this to float if your device does not support double computation
  22. typedef double value_type;
  23. //change this to host_vector< ... > of you want to run on CPU
  24. typedef thrust::device_vector< value_type > state_type;
  25. typedef thrust::device_vector< size_t > index_vector_type;
  26. // typedef thrust::host_vector< value_type > state_type;
  27. // typedef thrust::host_vector< size_t > index_vector_type;
  28. const value_type sigma = 10.0;
  29. const value_type b = 8.0 / 3.0;
  30. //[ thrust_lorenz_parameters_define_simple_system
  31. struct lorenz_system
  32. {
  33. struct lorenz_functor
  34. {
  35. template< class T >
  36. __host__ __device__
  37. void operator()( T t ) const
  38. {
  39. // unpack the parameter we want to vary and the Lorenz variables
  40. value_type R = thrust::get< 3 >( t );
  41. value_type x = thrust::get< 0 >( t );
  42. value_type y = thrust::get< 1 >( t );
  43. value_type z = thrust::get< 2 >( t );
  44. thrust::get< 4 >( t ) = sigma * ( y - x );
  45. thrust::get< 5 >( t ) = R * x - y - x * z;
  46. thrust::get< 6 >( t ) = -b * z + x * y ;
  47. }
  48. };
  49. lorenz_system( size_t N , const state_type &beta )
  50. : m_N( N ) , m_beta( beta ) { }
  51. template< class State , class Deriv >
  52. void operator()( const State &x , Deriv &dxdt , value_type t ) const
  53. {
  54. thrust::for_each(
  55. thrust::make_zip_iterator( thrust::make_tuple(
  56. boost::begin( x ) ,
  57. boost::begin( x ) + m_N ,
  58. boost::begin( x ) + 2 * m_N ,
  59. m_beta.begin() ,
  60. boost::begin( dxdt ) ,
  61. boost::begin( dxdt ) + m_N ,
  62. boost::begin( dxdt ) + 2 * m_N ) ) ,
  63. thrust::make_zip_iterator( thrust::make_tuple(
  64. boost::begin( x ) + m_N ,
  65. boost::begin( x ) + 2 * m_N ,
  66. boost::begin( x ) + 3 * m_N ,
  67. m_beta.begin() ,
  68. boost::begin( dxdt ) + m_N ,
  69. boost::begin( dxdt ) + 2 * m_N ,
  70. boost::begin( dxdt ) + 3 * m_N ) ) ,
  71. lorenz_functor() );
  72. }
  73. size_t m_N;
  74. const state_type &m_beta;
  75. };
  76. //]
  77. struct lorenz_perturbation_system
  78. {
  79. struct lorenz_perturbation_functor
  80. {
  81. template< class T >
  82. __host__ __device__
  83. void operator()( T t ) const
  84. {
  85. value_type R = thrust::get< 1 >( t );
  86. value_type x = thrust::get< 0 >( thrust::get< 0 >( t ) );
  87. value_type y = thrust::get< 1 >( thrust::get< 0 >( t ) );
  88. value_type z = thrust::get< 2 >( thrust::get< 0 >( t ) );
  89. value_type dx = thrust::get< 3 >( thrust::get< 0 >( t ) );
  90. value_type dy = thrust::get< 4 >( thrust::get< 0 >( t ) );
  91. value_type dz = thrust::get< 5 >( thrust::get< 0 >( t ) );
  92. thrust::get< 0 >( thrust::get< 2 >( t ) ) = sigma * ( y - x );
  93. thrust::get< 1 >( thrust::get< 2 >( t ) ) = R * x - y - x * z;
  94. thrust::get< 2 >( thrust::get< 2 >( t ) ) = -b * z + x * y ;
  95. thrust::get< 3 >( thrust::get< 2 >( t ) ) = sigma * ( dy - dx );
  96. thrust::get< 4 >( thrust::get< 2 >( t ) ) = ( R - z ) * dx - dy - x * dz;
  97. thrust::get< 5 >( thrust::get< 2 >( t ) ) = y * dx + x * dy - b * dz;
  98. }
  99. };
  100. lorenz_perturbation_system( size_t N , const state_type &beta )
  101. : m_N( N ) , m_beta( beta ) { }
  102. template< class State , class Deriv >
  103. void operator()( const State &x , Deriv &dxdt , value_type t ) const
  104. {
  105. thrust::for_each(
  106. thrust::make_zip_iterator( thrust::make_tuple(
  107. thrust::make_zip_iterator( thrust::make_tuple(
  108. boost::begin( x ) ,
  109. boost::begin( x ) + m_N ,
  110. boost::begin( x ) + 2 * m_N ,
  111. boost::begin( x ) + 3 * m_N ,
  112. boost::begin( x ) + 4 * m_N ,
  113. boost::begin( x ) + 5 * m_N ) ) ,
  114. m_beta.begin() ,
  115. thrust::make_zip_iterator( thrust::make_tuple(
  116. boost::begin( dxdt ) ,
  117. boost::begin( dxdt ) + m_N ,
  118. boost::begin( dxdt ) + 2 * m_N ,
  119. boost::begin( dxdt ) + 3 * m_N ,
  120. boost::begin( dxdt ) + 4 * m_N ,
  121. boost::begin( dxdt ) + 5 * m_N ) )
  122. ) ) ,
  123. thrust::make_zip_iterator( thrust::make_tuple(
  124. thrust::make_zip_iterator( thrust::make_tuple(
  125. boost::begin( x ) + m_N ,
  126. boost::begin( x ) + 2 * m_N ,
  127. boost::begin( x ) + 3 * m_N ,
  128. boost::begin( x ) + 4 * m_N ,
  129. boost::begin( x ) + 5 * m_N ,
  130. boost::begin( x ) + 6 * m_N ) ) ,
  131. m_beta.begin() ,
  132. thrust::make_zip_iterator( thrust::make_tuple(
  133. boost::begin( dxdt ) + m_N ,
  134. boost::begin( dxdt ) + 2 * m_N ,
  135. boost::begin( dxdt ) + 3 * m_N ,
  136. boost::begin( dxdt ) + 4 * m_N ,
  137. boost::begin( dxdt ) + 5 * m_N ,
  138. boost::begin( dxdt ) + 6 * m_N ) )
  139. ) ) ,
  140. lorenz_perturbation_functor() );
  141. }
  142. size_t m_N;
  143. const state_type &m_beta;
  144. };
  145. struct lyap_observer
  146. {
  147. //[thrust_lorenz_parameters_observer_functor
  148. struct lyap_functor
  149. {
  150. template< class T >
  151. __host__ __device__
  152. void operator()( T t ) const
  153. {
  154. value_type &dx = thrust::get< 0 >( t );
  155. value_type &dy = thrust::get< 1 >( t );
  156. value_type &dz = thrust::get< 2 >( t );
  157. value_type norm = sqrt( dx * dx + dy * dy + dz * dz );
  158. dx /= norm;
  159. dy /= norm;
  160. dz /= norm;
  161. thrust::get< 3 >( t ) += log( norm );
  162. }
  163. };
  164. //]
  165. lyap_observer( size_t N , size_t every = 100 )
  166. : m_N( N ) , m_lyap( N ) , m_every( every ) , m_count( 0 )
  167. {
  168. thrust::fill( m_lyap.begin() , m_lyap.end() , 0.0 );
  169. }
  170. template< class Lyap >
  171. void fill_lyap( Lyap &lyap )
  172. {
  173. thrust::copy( m_lyap.begin() , m_lyap.end() , lyap.begin() );
  174. for( size_t i=0 ; i<lyap.size() ; ++i )
  175. lyap[i] /= m_t_overall;
  176. }
  177. template< class State >
  178. void operator()( State &x , value_type t )
  179. {
  180. if( ( m_count != 0 ) && ( ( m_count % m_every ) == 0 ) )
  181. {
  182. thrust::for_each(
  183. thrust::make_zip_iterator( thrust::make_tuple(
  184. boost::begin( x ) + 3 * m_N ,
  185. boost::begin( x ) + 4 * m_N ,
  186. boost::begin( x ) + 5 * m_N ,
  187. m_lyap.begin() ) ) ,
  188. thrust::make_zip_iterator( thrust::make_tuple(
  189. boost::begin( x ) + 4 * m_N ,
  190. boost::begin( x ) + 5 * m_N ,
  191. boost::begin( x ) + 6 * m_N ,
  192. m_lyap.end() ) ) ,
  193. lyap_functor() );
  194. clog << t << "\n";
  195. }
  196. ++m_count;
  197. m_t_overall = t;
  198. }
  199. size_t m_N;
  200. state_type m_lyap;
  201. size_t m_every;
  202. size_t m_count;
  203. value_type m_t_overall;
  204. };
  205. const size_t N = 1024*2;
  206. const value_type dt = 0.01;
  207. int main( int arc , char* argv[] )
  208. {
  209. int driver_version , runtime_version;
  210. cudaDriverGetVersion( &driver_version );
  211. cudaRuntimeGetVersion ( &runtime_version );
  212. cout << driver_version << "\t" << runtime_version << endl;
  213. //[ thrust_lorenz_parameters_define_beta
  214. vector< value_type > beta_host( N );
  215. const value_type beta_min = 0.0 , beta_max = 56.0;
  216. for( size_t i=0 ; i<N ; ++i )
  217. beta_host[i] = beta_min + value_type( i ) * ( beta_max - beta_min ) / value_type( N - 1 );
  218. state_type beta = beta_host;
  219. //]
  220. //[ thrust_lorenz_parameters_integration
  221. state_type x( 6 * N );
  222. // initialize x,y,z
  223. thrust::fill( x.begin() , x.begin() + 3 * N , 10.0 );
  224. // initial dx
  225. thrust::fill( x.begin() + 3 * N , x.begin() + 4 * N , 1.0 );
  226. // initialize dy,dz
  227. thrust::fill( x.begin() + 4 * N , x.end() , 0.0 );
  228. // create error stepper, can be used with make_controlled or make_dense_output
  229. typedef runge_kutta_dopri5< state_type , value_type , state_type , value_type > stepper_type;
  230. lorenz_system lorenz( N , beta );
  231. lorenz_perturbation_system lorenz_perturbation( N , beta );
  232. lyap_observer obs( N , 1 );
  233. // calculate transients
  234. integrate_adaptive( make_controlled( 1.0e-6 , 1.0e-6 , stepper_type() ) , lorenz , std::make_pair( x.begin() , x.begin() + 3 * N ) , 0.0 , 10.0 , dt );
  235. // calculate the Lyapunov exponents -- the main loop
  236. double t = 0.0;
  237. while( t < 10000.0 )
  238. {
  239. integrate_adaptive( make_controlled( 1.0e-6 , 1.0e-6 , stepper_type() ) , lorenz_perturbation , x , t , t + 1.0 , 0.1 );
  240. t += 1.0;
  241. obs( x , t );
  242. }
  243. vector< value_type > lyap( N );
  244. obs.fill_lyap( lyap );
  245. for( size_t i=0 ; i<N ; ++i )
  246. cout << beta_host[i] << "\t" << lyap[i] << "\n";
  247. //]
  248. return 0;
  249. }