2 The example how the phase_oscillator ensemble can be implemented using CUDA and thrust
4 Copyright 2011-2013 Mario Mulansky
5 Copyright 2011 Karsten Ahnert
7 Distributed under the Boost Software License, Version 1.0.
8 (See accompanying file LICENSE_1_0.txt or
9 copy at http://www.boost.org/LICENSE_1_0.txt)
17 #include <thrust/device_vector.h>
18 #include <thrust/reduce.h>
19 #include <thrust/functional.h>
21 #include <boost/numeric/odeint.hpp>
22 #include <boost/numeric/odeint/external/thrust/thrust.hpp>
24 #include <boost/timer.hpp>
25 #include <boost/random/cauchy_distribution.hpp>
29 using namespace boost::numeric::odeint;
32 * Sorry for that dirty hack, but nvcc has large problems with boost::random.
34 * Nevertheless we need the cauchy distribution from boost::random, and therefore
35 * we need a generator. Here it is:
37 struct drand48_generator
39 typedef double result_type;
40 result_type operator()( void ) const { return drand48(); }
41 result_type min( void ) const { return 0.0; }
42 result_type max( void ) const { return 1.0; }
45 //[ thrust_phase_ensemble_state_type
46 //change this to float if your device does not support double computation
47 typedef double value_type;
49 //change this to host_vector< ... > of you want to run on CPU
50 typedef thrust::device_vector< value_type > state_type;
51 // typedef thrust::host_vector< value_type > state_type;
55 //[ thrust_phase_ensemble_mean_field_calculator
56 struct mean_field_calculator
58 struct sin_functor : public thrust::unary_function< value_type , value_type >
61 value_type operator()( value_type x) const
67 struct cos_functor : public thrust::unary_function< value_type , value_type >
70 value_type operator()( value_type x) const
76 static std::pair< value_type , value_type > get_mean( const state_type &x )
78 //[ thrust_phase_ensemble_sin_sum
79 value_type sin_sum = thrust::reduce(
80 thrust::make_transform_iterator( x.begin() , sin_functor() ) ,
81 thrust::make_transform_iterator( x.end() , sin_functor() ) );
83 value_type cos_sum = thrust::reduce(
84 thrust::make_transform_iterator( x.begin() , cos_functor() ) ,
85 thrust::make_transform_iterator( x.end() , cos_functor() ) );
87 cos_sum /= value_type( x.size() );
88 sin_sum /= value_type( x.size() );
90 value_type K = sqrt( cos_sum * cos_sum + sin_sum * sin_sum );
91 value_type Theta = atan2( sin_sum , cos_sum );
93 return std::make_pair( K , Theta );
100 //[ thrust_phase_ensemble_sys_function
101 class phase_oscillator_ensemble
108 value_type m_K , m_Theta , m_epsilon;
110 sys_functor( value_type K , value_type Theta , value_type epsilon )
111 : m_K( K ) , m_Theta( Theta ) , m_epsilon( epsilon ) { }
113 template< class Tuple >
115 void operator()( Tuple t )
117 thrust::get<2>(t) = thrust::get<1>(t) + m_epsilon * m_K * sin( m_Theta - thrust::get<0>(t) );
123 phase_oscillator_ensemble( size_t N , value_type g = 1.0 , value_type epsilon = 1.0 )
124 : m_omega() , m_N( N ) , m_epsilon( epsilon )
126 create_frequencies( g );
129 void create_frequencies( value_type g )
131 boost::cauchy_distribution< value_type > cauchy( 0.0 , g );
132 // boost::variate_generator< boost::mt19937&, boost::cauchy_distribution< value_type > > gen( rng , cauchy );
133 drand48_generator d48;
134 vector< value_type > omega( m_N );
135 for( size_t i=0 ; i<m_N ; ++i )
136 omega[i] = cauchy( d48 );
137 // generate( omega.begin() , omega.end() , gen );
141 void set_epsilon( value_type epsilon ) { m_epsilon = epsilon; }
143 value_type get_epsilon( void ) const { return m_epsilon; }
146 void operator() ( const state_type &x , state_type &dxdt , const value_type dt ) const
148 std::pair< value_type , value_type > mean_field = mean_field_calculator::get_mean( x );
151 thrust::make_zip_iterator( thrust::make_tuple( x.begin() , m_omega.begin() , dxdt.begin() ) ),
152 thrust::make_zip_iterator( thrust::make_tuple( x.end() , m_omega.end() , dxdt.end()) ) ,
153 sys_functor( mean_field.first , mean_field.second , m_epsilon )
163 value_type m_epsilon;
169 //[ thrust_phase_ensemble_observer
170 struct statistics_observer
175 statistics_observer( void )
176 : m_K_mean( 0.0 ) , m_count( 0 ) { }
178 template< class State >
179 void operator()( const State &x , value_type t )
181 std::pair< value_type , value_type > mean = mean_field_calculator::get_mean( x );
182 m_K_mean += mean.first;
186 value_type get_K_mean( void ) const { return ( m_count != 0 ) ? m_K_mean / value_type( m_count ) : 0.0 ; }
188 void reset( void ) { m_K_mean = 0.0; m_count = 0; }
194 // const size_t N = 16384 * 128;
195 const size_t N = 16384;
196 const value_type pi = 3.1415926535897932384626433832795029;
197 const value_type dt = 0.1;
198 const value_type d_epsilon = 0.1;
199 const value_type epsilon_min = 0.0;
200 const value_type epsilon_max = 5.0;
201 const value_type t_transients = 10.0;
202 const value_type t_max = 100.0;
204 int main( int arc , char* argv[] )
206 // initial conditions on host
207 vector< value_type > x_host( N );
208 for( size_t i=0 ; i<N ; ++i ) x_host[i] = 2.0 * pi * drand48();
210 //[ thrust_phase_ensemble_system_instance
211 phase_oscillator_ensemble ensemble( N , 1.0 );
217 boost::timer timer_local;
218 double dopri5_time = 0.0 , rk4_time = 0.0;
220 //[thrust_phase_ensemble_define_dopri5
221 typedef runge_kutta_dopri5< state_type , value_type , state_type , value_type > stepper_type;
224 ofstream fout( "phase_ensemble_dopri5.dat" );
226 for( value_type epsilon = epsilon_min ; epsilon < epsilon_max ; epsilon += d_epsilon )
228 ensemble.set_epsilon( epsilon );
229 statistics_observer obs;
230 state_type x = x_host;
232 timer_local.restart();
234 // calculate some transients steps
235 //[ thrust_phase_ensemble_integration
236 size_t steps1 = integrate_const( make_controlled( 1.0e-6 , 1.0e-6 , stepper_type() ) , boost::ref( ensemble ) , x , 0.0 , t_transients , dt );
239 // integrate and compute the statistics
240 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 ) );
242 fout << epsilon << "\t" << obs.get_K_mean() << endl;
243 cout << "Dopri5 : " << epsilon << "\t" << obs.get_K_mean() << "\t" << timer_local.elapsed() << "\t" << steps1 << "\t" << steps2 << endl;
245 dopri5_time = timer.elapsed();
251 //[ thrust_phase_ensemble_define_rk4
252 typedef runge_kutta4< state_type , value_type , state_type , value_type > stepper_type;
255 ofstream fout( "phase_ensemble_rk4.dat" );
257 for( value_type epsilon = epsilon_min ; epsilon < epsilon_max ; epsilon += d_epsilon )
259 ensemble.set_epsilon( epsilon );
260 statistics_observer obs;
261 state_type x = x_host;
263 timer_local.restart();
265 // calculate some transients steps
266 size_t steps1 = integrate_const( stepper_type() , boost::ref( ensemble ) , x , 0.0 , t_transients , dt );
268 // integrate and compute the statistics
269 size_t steps2 = integrate_const( stepper_type() , boost::ref( ensemble ) , x , 0.0 , t_max , dt , boost::ref( obs ) );
270 fout << epsilon << "\t" << obs.get_K_mean() << endl;
271 cout << "RK4 : " << epsilon << "\t" << obs.get_K_mean() << "\t" << timer_local.elapsed() << "\t" << steps1 << "\t" << steps2 << endl;
273 rk4_time = timer.elapsed();
276 cout << "Dopri 5 : " << dopri5_time << " s\n";
277 cout << "RK4 : " << rk4_time << "\n";