• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
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 
13 #include <iostream>
14 #include <map>
15 #include <string>
16 #include <fstream>
17 
18 #ifndef M_PI //not there on windows
19 #define M_PI 3.1415927 //...
20 #endif
21 
22 #include <boost/numeric/odeint.hpp>
23 
24 using namespace std;
25 using namespace boost::numeric::odeint;
26 
27 //[ two_dimensional_phase_lattice_definition
28 typedef boost::numeric::ublas::matrix< double > state_type;
29 
30 struct two_dimensional_phase_lattice
31 {
two_dimensional_phase_latticetwo_dimensional_phase_lattice32     two_dimensional_phase_lattice( double gamma = 0.5 )
33     : m_gamma( gamma ) { }
34 
operator ()two_dimensional_phase_lattice35     void operator()( const state_type &x , state_type &dxdt , double /* t */ ) const
36     {
37         size_t size1 = x.size1() , size2 = x.size2();
38 
39         for( size_t i=1 ; i<size1-1 ; ++i )
40         {
41             for( size_t j=1 ; j<size2-1 ; ++j )
42             {
43                 dxdt( i , j ) =
44                         coupling_func( x( i + 1 , j ) - x( i , j ) ) +
45                         coupling_func( x( i - 1 , j ) - x( i , j ) ) +
46                         coupling_func( x( i , j + 1 ) - x( i , j ) ) +
47                         coupling_func( x( i , j - 1 ) - x( i , j ) );
48             }
49         }
50 
51         for( size_t i=0 ; i<x.size1() ; ++i ) dxdt( i , 0 ) = dxdt( i , x.size2() -1 ) = 0.0;
52         for( size_t j=0 ; j<x.size2() ; ++j ) dxdt( 0 , j ) = dxdt( x.size1() -1 , j ) = 0.0;
53     }
54 
coupling_functwo_dimensional_phase_lattice55     double coupling_func( double x ) const
56     {
57         return sin( x ) - m_gamma * ( 1.0 - cos( x ) );
58     }
59 
60     double m_gamma;
61 };
62 //]
63 
64 
65 struct write_for_gnuplot
66 {
67     size_t m_every , m_count;
68 
write_for_gnuplotwrite_for_gnuplot69     write_for_gnuplot( size_t every = 10 )
70     : m_every( every ) , m_count( 0 ) { }
71 
operator ()write_for_gnuplot72     void operator()( const state_type &x , double t )
73     {
74         if( ( m_count % m_every ) == 0 )
75         {
76             clog << t << endl;
77             cout << "sp '-'" << endl;
78             for( size_t i=0 ; i<x.size1() ; ++i )
79             {
80                 for( size_t j=0 ; j<x.size2() ; ++j )
81                 {
82                     cout << i << "\t" << j << "\t" << sin( x( i , j ) ) << "\n";
83                 }
84                 cout << "\n";
85             }
86             cout << "e" << endl;
87         }
88 
89         ++m_count;
90     }
91 };
92 
93 class write_snapshots
94 {
95 public:
96 
97     typedef std::map< size_t , std::string > map_type;
98 
write_snapshots(void)99     write_snapshots( void ) : m_count( 0 ) { }
100 
operator ()(const state_type & x,double t)101     void operator()( const state_type &x , double t )
102     {
103         map< size_t , string >::const_iterator it = m_snapshots.find( m_count );
104         if( it != m_snapshots.end() )
105         {
106             ofstream fout( it->second.c_str() );
107             for( size_t i=0 ; i<x.size1() ; ++i )
108             {
109                 for( size_t j=0 ; j<x.size2() ; ++j )
110                 {
111                     fout << i << "\t" << j << "\t" << x( i , j ) << "\t" << sin( x( i , j ) ) << "\n";
112                 }
113                 fout << "\n";
114             }
115         }
116         ++m_count;
117     }
118 
snapshots(void)119     map_type& snapshots( void ) { return m_snapshots; }
snapshots(void) const120     const map_type& snapshots( void ) const { return m_snapshots; }
121 
122 private:
123 
124     size_t m_count;
125     map_type m_snapshots;
126 };
127 
128 
main(int argc,char ** argv)129 int main( int argc , char **argv )
130 {
131     size_t size1 = 128 , size2 = 128;
132     state_type x( size1 , size2 , 0.0 );
133 
134     for( size_t i=(size1/2-10) ; i<(size1/2+10) ; ++i )
135         for( size_t j=(size2/2-10) ; j<(size2/2+10) ; ++j )
136             x( i , j ) = static_cast<double>( rand() ) / RAND_MAX * 2.0 * M_PI;
137 
138     write_snapshots snapshots;
139     snapshots.snapshots().insert( make_pair( size_t( 0 ) , string( "lat_0000.dat" ) ) );
140     snapshots.snapshots().insert( make_pair( size_t( 100 ) , string( "lat_0100.dat" ) ) );
141     snapshots.snapshots().insert( make_pair( size_t( 1000 ) , string( "lat_1000.dat" ) ) );
142     observer_collection< state_type , double > obs;
143     obs.observers().push_back( write_for_gnuplot( 10 ) );
144     obs.observers().push_back( snapshots );
145 
146     cout << "set term x11" << endl;
147     cout << "set pm3d map" << endl;
148 
149     integrate_const( runge_kutta4<state_type>() , two_dimensional_phase_lattice( 1.2 ) ,
150                      x , 0.0 , 1001.0 , 0.1 , boost::ref( obs ) );
151 
152     // controlled steppers work only after ublas bugfix
153     //integrate_const( make_dense_output< runge_kutta_dopri5< state_type > >( 1E-6 , 1E-6 ) , two_dimensional_phase_lattice( 1.2 ) ,
154     //        x , 0.0 , 1001.0 , 0.1 , boost::ref( obs ) );
155 
156 
157     return 0;
158 }
159