phase_oscillator_ensemble.cu 8.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280
  1. /*
  2. The example how the phase_oscillator ensemble can be implemented using CUDA and thrust
  3. Copyright 2011-2013 Mario Mulansky
  4. Copyright 2011 Karsten Ahnert
  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. #include <iostream>
  10. #include <fstream>
  11. #include <cmath>
  12. #include <utility>
  13. #include <thrust/device_vector.h>
  14. #include <thrust/reduce.h>
  15. #include <thrust/functional.h>
  16. #include <boost/numeric/odeint.hpp>
  17. #include <boost/numeric/odeint/external/thrust/thrust.hpp>
  18. #include <boost/timer.hpp>
  19. #include <boost/random/cauchy_distribution.hpp>
  20. using namespace std;
  21. using namespace boost::numeric::odeint;
  22. /*
  23. * Sorry for that dirty hack, but nvcc has large problems with boost::random.
  24. *
  25. * Nevertheless we need the cauchy distribution from boost::random, and therefore
  26. * we need a generator. Here it is:
  27. */
  28. struct drand48_generator
  29. {
  30. typedef double result_type;
  31. result_type operator()( void ) const { return drand48(); }
  32. result_type min( void ) const { return 0.0; }
  33. result_type max( void ) const { return 1.0; }
  34. };
  35. //[ thrust_phase_ensemble_state_type
  36. //change this to float if your device does not support double computation
  37. typedef double value_type;
  38. //change this to host_vector< ... > of you want to run on CPU
  39. typedef thrust::device_vector< value_type > state_type;
  40. // typedef thrust::host_vector< value_type > state_type;
  41. //]
  42. //[ thrust_phase_ensemble_mean_field_calculator
  43. struct mean_field_calculator
  44. {
  45. struct sin_functor : public thrust::unary_function< value_type , value_type >
  46. {
  47. __host__ __device__
  48. value_type operator()( value_type x) const
  49. {
  50. return sin( x );
  51. }
  52. };
  53. struct cos_functor : public thrust::unary_function< value_type , value_type >
  54. {
  55. __host__ __device__
  56. value_type operator()( value_type x) const
  57. {
  58. return cos( x );
  59. }
  60. };
  61. static std::pair< value_type , value_type > get_mean( const state_type &x )
  62. {
  63. //[ thrust_phase_ensemble_sin_sum
  64. value_type sin_sum = thrust::reduce(
  65. thrust::make_transform_iterator( x.begin() , sin_functor() ) ,
  66. thrust::make_transform_iterator( x.end() , sin_functor() ) );
  67. //]
  68. value_type cos_sum = thrust::reduce(
  69. thrust::make_transform_iterator( x.begin() , cos_functor() ) ,
  70. thrust::make_transform_iterator( x.end() , cos_functor() ) );
  71. cos_sum /= value_type( x.size() );
  72. sin_sum /= value_type( x.size() );
  73. value_type K = sqrt( cos_sum * cos_sum + sin_sum * sin_sum );
  74. value_type Theta = atan2( sin_sum , cos_sum );
  75. return std::make_pair( K , Theta );
  76. }
  77. };
  78. //]
  79. //[ thrust_phase_ensemble_sys_function
  80. class phase_oscillator_ensemble
  81. {
  82. public:
  83. struct sys_functor
  84. {
  85. value_type m_K , m_Theta , m_epsilon;
  86. sys_functor( value_type K , value_type Theta , value_type epsilon )
  87. : m_K( K ) , m_Theta( Theta ) , m_epsilon( epsilon ) { }
  88. template< class Tuple >
  89. __host__ __device__
  90. void operator()( Tuple t )
  91. {
  92. thrust::get<2>(t) = thrust::get<1>(t) + m_epsilon * m_K * sin( m_Theta - thrust::get<0>(t) );
  93. }
  94. };
  95. // ...
  96. //<-
  97. phase_oscillator_ensemble( size_t N , value_type g = 1.0 , value_type epsilon = 1.0 )
  98. : m_omega() , m_N( N ) , m_epsilon( epsilon )
  99. {
  100. create_frequencies( g );
  101. }
  102. void create_frequencies( value_type g )
  103. {
  104. boost::cauchy_distribution< value_type > cauchy( 0.0 , g );
  105. // boost::variate_generator< boost::mt19937&, boost::cauchy_distribution< value_type > > gen( rng , cauchy );
  106. drand48_generator d48;
  107. vector< value_type > omega( m_N );
  108. for( size_t i=0 ; i<m_N ; ++i )
  109. omega[i] = cauchy( d48 );
  110. // generate( omega.begin() , omega.end() , gen );
  111. m_omega = omega;
  112. }
  113. void set_epsilon( value_type epsilon ) { m_epsilon = epsilon; }
  114. value_type get_epsilon( void ) const { return m_epsilon; }
  115. //->
  116. void operator() ( const state_type &x , state_type &dxdt , const value_type dt ) const
  117. {
  118. std::pair< value_type , value_type > mean_field = mean_field_calculator::get_mean( x );
  119. thrust::for_each(
  120. thrust::make_zip_iterator( thrust::make_tuple( x.begin() , m_omega.begin() , dxdt.begin() ) ),
  121. thrust::make_zip_iterator( thrust::make_tuple( x.end() , m_omega.end() , dxdt.end()) ) ,
  122. sys_functor( mean_field.first , mean_field.second , m_epsilon )
  123. );
  124. }
  125. // ...
  126. //<-
  127. private:
  128. state_type m_omega;
  129. const size_t m_N;
  130. value_type m_epsilon;
  131. //->
  132. };
  133. //]
  134. //[ thrust_phase_ensemble_observer
  135. struct statistics_observer
  136. {
  137. value_type m_K_mean;
  138. size_t m_count;
  139. statistics_observer( void )
  140. : m_K_mean( 0.0 ) , m_count( 0 ) { }
  141. template< class State >
  142. void operator()( const State &x , value_type t )
  143. {
  144. std::pair< value_type , value_type > mean = mean_field_calculator::get_mean( x );
  145. m_K_mean += mean.first;
  146. ++m_count;
  147. }
  148. value_type get_K_mean( void ) const { return ( m_count != 0 ) ? m_K_mean / value_type( m_count ) : 0.0 ; }
  149. void reset( void ) { m_K_mean = 0.0; m_count = 0; }
  150. };
  151. //]
  152. // const size_t N = 16384 * 128;
  153. const size_t N = 16384;
  154. const value_type pi = 3.1415926535897932384626433832795029;
  155. const value_type dt = 0.1;
  156. const value_type d_epsilon = 0.1;
  157. const value_type epsilon_min = 0.0;
  158. const value_type epsilon_max = 5.0;
  159. const value_type t_transients = 10.0;
  160. const value_type t_max = 100.0;
  161. int main( int arc , char* argv[] )
  162. {
  163. // initial conditions on host
  164. vector< value_type > x_host( N );
  165. for( size_t i=0 ; i<N ; ++i ) x_host[i] = 2.0 * pi * drand48();
  166. //[ thrust_phase_ensemble_system_instance
  167. phase_oscillator_ensemble ensemble( N , 1.0 );
  168. //]
  169. boost::timer timer;
  170. boost::timer timer_local;
  171. double dopri5_time = 0.0 , rk4_time = 0.0;
  172. {
  173. //[thrust_phase_ensemble_define_dopri5
  174. typedef runge_kutta_dopri5< state_type , value_type , state_type , value_type > stepper_type;
  175. //]
  176. ofstream fout( "phase_ensemble_dopri5.dat" );
  177. timer.restart();
  178. for( value_type epsilon = epsilon_min ; epsilon < epsilon_max ; epsilon += d_epsilon )
  179. {
  180. ensemble.set_epsilon( epsilon );
  181. statistics_observer obs;
  182. state_type x = x_host;
  183. timer_local.restart();
  184. // calculate some transients steps
  185. //[ thrust_phase_ensemble_integration
  186. size_t steps1 = integrate_const( make_controlled( 1.0e-6 , 1.0e-6 , stepper_type() ) , boost::ref( ensemble ) , x , 0.0 , t_transients , dt );
  187. //]
  188. // integrate and compute the statistics
  189. size_t steps2 = integrate_const( make_dense_output( 1.0e-6 , 1.0e-6 , stepper_type() ) , boost::ref( ensemble ) , x , 0.0 , t_max , dt , boost::ref( obs ) );
  190. fout << epsilon << "\t" << obs.get_K_mean() << endl;
  191. cout << "Dopri5 : " << epsilon << "\t" << obs.get_K_mean() << "\t" << timer_local.elapsed() << "\t" << steps1 << "\t" << steps2 << endl;
  192. }
  193. dopri5_time = timer.elapsed();
  194. }
  195. {
  196. //[ thrust_phase_ensemble_define_rk4
  197. typedef runge_kutta4< state_type , value_type , state_type , value_type > stepper_type;
  198. //]
  199. ofstream fout( "phase_ensemble_rk4.dat" );
  200. timer.restart();
  201. for( value_type epsilon = epsilon_min ; epsilon < epsilon_max ; epsilon += d_epsilon )
  202. {
  203. ensemble.set_epsilon( epsilon );
  204. statistics_observer obs;
  205. state_type x = x_host;
  206. timer_local.restart();
  207. // calculate some transients steps
  208. size_t steps1 = integrate_const( stepper_type() , boost::ref( ensemble ) , x , 0.0 , t_transients , dt );
  209. // integrate and compute the statistics
  210. size_t steps2 = integrate_const( stepper_type() , boost::ref( ensemble ) , x , 0.0 , t_max , dt , boost::ref( obs ) );
  211. fout << epsilon << "\t" << obs.get_K_mean() << endl;
  212. cout << "RK4 : " << epsilon << "\t" << obs.get_K_mean() << "\t" << timer_local.elapsed() << "\t" << steps1 << "\t" << steps2 << endl;
  213. }
  214. rk4_time = timer.elapsed();
  215. }
  216. cout << "Dopri 5 : " << dopri5_time << " s\n";
  217. cout << "RK4 : " << rk4_time << "\n";
  218. return 0;
  219. }