phase_chain.cpp 2.5 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495
  1. /*
  2. * phase_chain.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. //[phase_chain_openmp_header
  18. #include <omp.h>
  19. #include <boost/numeric/odeint.hpp>
  20. #include <boost/numeric/odeint/external/openmp/openmp.hpp>
  21. //]
  22. using namespace std;
  23. using namespace boost::numeric::odeint;
  24. using boost::timer::cpu_timer;
  25. using boost::math::double_constants::pi;
  26. //[phase_chain_vector_state
  27. typedef std::vector< double > state_type;
  28. //]
  29. //[phase_chain_rhs
  30. struct phase_chain
  31. {
  32. phase_chain( double gamma = 0.5 )
  33. : m_gamma( gamma ) { }
  34. void operator()( const state_type &x , state_type &dxdt , double /* t */ ) const
  35. {
  36. const size_t N = x.size();
  37. #pragma omp parallel for schedule(runtime)
  38. for(size_t i = 1 ; i < N - 1 ; ++i)
  39. {
  40. dxdt[i] = coupling_func( x[i+1] - x[i] ) +
  41. coupling_func( x[i-1] - x[i] );
  42. }
  43. dxdt[0 ] = coupling_func( x[1 ] - x[0 ] );
  44. dxdt[N-1] = coupling_func( x[N-2] - x[N-1] );
  45. }
  46. double coupling_func( double x ) const
  47. {
  48. return sin( x ) - m_gamma * ( 1.0 - cos( x ) );
  49. }
  50. double m_gamma;
  51. };
  52. //]
  53. int main( int argc , char **argv )
  54. {
  55. //[phase_chain_init
  56. size_t N = 131101;
  57. state_type x( N );
  58. boost::random::uniform_real_distribution<double> distribution( 0.0 , 2.0*pi );
  59. boost::random::mt19937 engine( 0 );
  60. generate( x.begin() , x.end() , boost::bind( distribution , engine ) );
  61. //]
  62. //[phase_chain_stepper
  63. typedef runge_kutta4<
  64. state_type , double ,
  65. state_type , double ,
  66. openmp_range_algebra
  67. > stepper_type;
  68. //]
  69. //[phase_chain_scheduling
  70. int chunk_size = N/omp_get_max_threads();
  71. omp_set_schedule( omp_sched_static , chunk_size );
  72. //]
  73. cpu_timer timer;
  74. //[phase_chain_integrate
  75. integrate_n_steps( stepper_type() , phase_chain( 1.2 ) ,
  76. x , 0.0 , 0.01 , 100 );
  77. //]
  78. double run_time = static_cast<double>(timer.elapsed().wall) * 1.0e-9;
  79. std::cerr << run_time << "s" << std::endl;
  80. // copy(x.begin(), x.end(), ostream_iterator<double>(cout, "\n"));
  81. return 0;
  82. }