two_dimensional_phase_lattice.cpp 4.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158
  1. /*
  2. * two_dimensional_phase_lattice.cpp
  3. *
  4. * This example show how one can use matrices as state types in odeint.
  5. *
  6. * Copyright 2011-2012 Karsten Ahnert
  7. * Copyright 2011-2013 Mario Mulansky
  8. * Distributed under the Boost Software License, Version 1.0. (See
  9. * accompanying file LICENSE_1_0.txt or copy at
  10. * http://www.boost.org/LICENSE_1_0.txt)
  11. */
  12. #include <iostream>
  13. #include <map>
  14. #include <string>
  15. #include <fstream>
  16. #ifndef M_PI //not there on windows
  17. #define M_PI 3.1415927 //...
  18. #endif
  19. #include <boost/numeric/odeint.hpp>
  20. using namespace std;
  21. using namespace boost::numeric::odeint;
  22. //[ two_dimensional_phase_lattice_definition
  23. typedef boost::numeric::ublas::matrix< double > state_type;
  24. struct two_dimensional_phase_lattice
  25. {
  26. two_dimensional_phase_lattice( double gamma = 0.5 )
  27. : m_gamma( gamma ) { }
  28. void operator()( const state_type &x , state_type &dxdt , double /* t */ ) const
  29. {
  30. size_t size1 = x.size1() , size2 = x.size2();
  31. for( size_t i=1 ; i<size1-1 ; ++i )
  32. {
  33. for( size_t j=1 ; j<size2-1 ; ++j )
  34. {
  35. dxdt( i , j ) =
  36. coupling_func( x( i + 1 , j ) - x( i , j ) ) +
  37. coupling_func( x( i - 1 , j ) - x( i , j ) ) +
  38. coupling_func( x( i , j + 1 ) - x( i , j ) ) +
  39. coupling_func( x( i , j - 1 ) - x( i , j ) );
  40. }
  41. }
  42. for( size_t i=0 ; i<x.size1() ; ++i ) dxdt( i , 0 ) = dxdt( i , x.size2() -1 ) = 0.0;
  43. for( size_t j=0 ; j<x.size2() ; ++j ) dxdt( 0 , j ) = dxdt( x.size1() -1 , j ) = 0.0;
  44. }
  45. double coupling_func( double x ) const
  46. {
  47. return sin( x ) - m_gamma * ( 1.0 - cos( x ) );
  48. }
  49. double m_gamma;
  50. };
  51. //]
  52. struct write_for_gnuplot
  53. {
  54. size_t m_every , m_count;
  55. write_for_gnuplot( size_t every = 10 )
  56. : m_every( every ) , m_count( 0 ) { }
  57. void operator()( const state_type &x , double t )
  58. {
  59. if( ( m_count % m_every ) == 0 )
  60. {
  61. clog << t << endl;
  62. cout << "sp '-'" << endl;
  63. for( size_t i=0 ; i<x.size1() ; ++i )
  64. {
  65. for( size_t j=0 ; j<x.size2() ; ++j )
  66. {
  67. cout << i << "\t" << j << "\t" << sin( x( i , j ) ) << "\n";
  68. }
  69. cout << "\n";
  70. }
  71. cout << "e" << endl;
  72. }
  73. ++m_count;
  74. }
  75. };
  76. class write_snapshots
  77. {
  78. public:
  79. typedef std::map< size_t , std::string > map_type;
  80. write_snapshots( void ) : m_count( 0 ) { }
  81. void operator()( const state_type &x , double t )
  82. {
  83. map< size_t , string >::const_iterator it = m_snapshots.find( m_count );
  84. if( it != m_snapshots.end() )
  85. {
  86. ofstream fout( it->second.c_str() );
  87. for( size_t i=0 ; i<x.size1() ; ++i )
  88. {
  89. for( size_t j=0 ; j<x.size2() ; ++j )
  90. {
  91. fout << i << "\t" << j << "\t" << x( i , j ) << "\t" << sin( x( i , j ) ) << "\n";
  92. }
  93. fout << "\n";
  94. }
  95. }
  96. ++m_count;
  97. }
  98. map_type& snapshots( void ) { return m_snapshots; }
  99. const map_type& snapshots( void ) const { return m_snapshots; }
  100. private:
  101. size_t m_count;
  102. map_type m_snapshots;
  103. };
  104. int main( int argc , char **argv )
  105. {
  106. size_t size1 = 128 , size2 = 128;
  107. state_type x( size1 , size2 , 0.0 );
  108. for( size_t i=(size1/2-10) ; i<(size1/2+10) ; ++i )
  109. for( size_t j=(size2/2-10) ; j<(size2/2+10) ; ++j )
  110. x( i , j ) = static_cast<double>( rand() ) / RAND_MAX * 2.0 * M_PI;
  111. write_snapshots snapshots;
  112. snapshots.snapshots().insert( make_pair( size_t( 0 ) , string( "lat_0000.dat" ) ) );
  113. snapshots.snapshots().insert( make_pair( size_t( 100 ) , string( "lat_0100.dat" ) ) );
  114. snapshots.snapshots().insert( make_pair( size_t( 1000 ) , string( "lat_1000.dat" ) ) );
  115. observer_collection< state_type , double > obs;
  116. obs.observers().push_back( write_for_gnuplot( 10 ) );
  117. obs.observers().push_back( snapshots );
  118. cout << "set term x11" << endl;
  119. cout << "set pm3d map" << endl;
  120. integrate_const( runge_kutta4<state_type>() , two_dimensional_phase_lattice( 1.2 ) ,
  121. x , 0.0 , 1001.0 , 0.1 , boost::ref( obs ) );
  122. // controlled steppers work only after ublas bugfix
  123. //integrate_const( make_dense_output< runge_kutta_dopri5< state_type > >( 1E-6 , 1E-6 ) , two_dimensional_phase_lattice( 1.2 ) ,
  124. // x , 0.0 , 1001.0 , 0.1 , boost::ref( obs ) );
  125. return 0;
  126. }