phase_chain_omp_state.cpp 2.7 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798
  1. /*
  2. * phase_chain_omp_state.cpp
  3. *
  4. * Example of OMP parallelization with odeint
  5. *
  6. * Copyright 2013 Karsten Ahnert
  7. * Copyright 2013 Mario Mulansky
  8. * Copyright 2013 Pascal Germroth
  9. * Distributed under the Boost Software License, Version 1.0. (See
  10. * accompanying file LICENSE_1_0.txt or copy at
  11. * http://www.boost.org/LICENSE_1_0.txt)
  12. */
  13. #include <iostream>
  14. #include <vector>
  15. #include <boost/random.hpp>
  16. #include <boost/timer/timer.hpp>
  17. #include <omp.h>
  18. #include <boost/numeric/odeint.hpp>
  19. #include <boost/numeric/odeint/external/openmp/openmp.hpp>
  20. #include <boost/numeric/odeint.hpp>
  21. using namespace std;
  22. using namespace boost::numeric::odeint;
  23. using boost::timer::cpu_timer;
  24. using boost::math::double_constants::pi;
  25. typedef openmp_state<double> state_type;
  26. //[phase_chain_state_rhs
  27. struct phase_chain_omp_state
  28. {
  29. phase_chain_omp_state( double gamma = 0.5 )
  30. : m_gamma( gamma ) { }
  31. void operator()( const state_type &x , state_type &dxdt , double /* t */ ) const
  32. {
  33. const size_t N = x.size();
  34. #pragma omp parallel for schedule(runtime)
  35. for(size_t n = 0 ; n < N ; ++n)
  36. {
  37. const size_t M = x[n].size();
  38. for(size_t m = 1 ; m < M-1 ; ++m)
  39. {
  40. dxdt[n][m] = coupling_func( x[n][m+1] - x[n][m] ) +
  41. coupling_func( x[n][m-1] - x[n][m] );
  42. }
  43. dxdt[n][0] = coupling_func( x[n][1] - x[n][0] );
  44. if( n > 0 )
  45. {
  46. dxdt[n][0] += coupling_func( x[n-1].back() - x[n].front() );
  47. }
  48. dxdt[n][M-1] = coupling_func( x[n][M-2] - x[n][M-1] );
  49. if( n < N-1 )
  50. {
  51. dxdt[n][M-1] += coupling_func( x[n+1].front() - x[n].back() );
  52. }
  53. }
  54. }
  55. double coupling_func( double x ) const
  56. {
  57. return sin( x ) - m_gamma * ( 1.0 - cos( x ) );
  58. }
  59. double m_gamma;
  60. };
  61. //]
  62. int main( int argc , char **argv )
  63. {
  64. //[phase_chain_state_init
  65. const size_t N = 131101;
  66. vector<double> x( N );
  67. boost::random::uniform_real_distribution<double> distribution( 0.0 , 2.0*pi );
  68. boost::random::mt19937 engine( 0 );
  69. generate( x.begin() , x.end() , boost::bind( distribution , engine ) );
  70. const size_t blocks = omp_get_max_threads();
  71. state_type x_split( blocks );
  72. split( x , x_split );
  73. //]
  74. cpu_timer timer;
  75. //[phase_chain_state_integrate
  76. integrate_n_steps( runge_kutta4<state_type>() , phase_chain_omp_state( 1.2 ) ,
  77. x_split , 0.0 , 0.01 , 100 );
  78. unsplit( x_split , x );
  79. //]
  80. double run_time = static_cast<double>(timer.elapsed().wall) * 1.0e-9;
  81. std::cerr << run_time << "s" << std::endl;
  82. // copy(x.begin(), x.end(), ostream_iterator<double>(cout, "\n"));
  83. return 0;
  84. }