relaxation.cu 2.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081
  1. /*
  2. Copyright 2011-2012 Karsten Ahnert
  3. Copyright 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. /*
  9. * Solves many relaxation equations dxdt = - a * x in parallel and for different values of a.
  10. * The relaxation equations are completely uncoupled.
  11. */
  12. #include <thrust/device_vector.h>
  13. #include <boost/ref.hpp>
  14. #include <boost/numeric/odeint.hpp>
  15. #include <boost/numeric/odeint/external/thrust/thrust.hpp>
  16. using namespace std;
  17. using namespace boost::numeric::odeint;
  18. // change to float if your GPU does not support doubles
  19. typedef double value_type;
  20. typedef thrust::device_vector< value_type > state_type;
  21. typedef runge_kutta4< state_type , value_type , state_type , value_type > stepper_type;
  22. struct relaxation
  23. {
  24. struct relaxation_functor
  25. {
  26. template< class T >
  27. __host__ __device__
  28. void operator()( T t ) const
  29. {
  30. // unpack the parameter we want to vary and the Lorenz variables
  31. value_type a = thrust::get< 1 >( t );
  32. value_type x = thrust::get< 0 >( t );
  33. thrust::get< 2 >( t ) = -a * x;
  34. }
  35. };
  36. relaxation( size_t N , const state_type &a )
  37. : m_N( N ) , m_a( a ) { }
  38. void operator()( const state_type &x , state_type &dxdt , value_type t ) const
  39. {
  40. thrust::for_each(
  41. thrust::make_zip_iterator( thrust::make_tuple( x.begin() , m_a.begin() , dxdt.begin() ) ) ,
  42. thrust::make_zip_iterator( thrust::make_tuple( x.end() , m_a.end() , dxdt.end() ) ) ,
  43. relaxation_functor() );
  44. }
  45. size_t m_N;
  46. const state_type &m_a;
  47. };
  48. const size_t N = 1024 * 1024;
  49. const value_type dt = 0.01;
  50. int main( int arc , char* argv[] )
  51. {
  52. // initialize the relaxation constants a
  53. vector< value_type > a_host( N );
  54. for( size_t i=0 ; i<N ; ++i ) a_host[i] = drand48();
  55. state_type a = a_host;
  56. // initialize the intial state x
  57. state_type x( N );
  58. thrust::fill( x.begin() , x.end() , 1.0 );
  59. // integrate
  60. relaxation relax( N , a );
  61. integrate_const( stepper_type() , boost::ref( relax ) , x , 0.0 , 10.0 , dt );
  62. return 0;
  63. }