-
Notifications
You must be signed in to change notification settings - Fork 102
/
harmonic_oscillator.cpp
211 lines (136 loc) · 4.58 KB
/
harmonic_oscillator.cpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
/*
Copyright 2010-2012 Karsten Ahnert
Copyright 2011-2013 Mario Mulansky
Copyright 2013 Pascal Germroth
Distributed under the Boost Software License, Version 1.0.
(See accompanying file LICENSE_1_0.txt or
copy at http://www.boost.org/LICENSE_1_0.txt)
*/
#include <iostream>
#include <vector>
#include <boost/numeric/odeint.hpp>
//[ rhs_function
/* The type of container used to hold the state vector */
typedef std::vector< double > state_type;
const double gam = 0.15;
/* The rhs of x' = f(x) */
void harmonic_oscillator( const state_type &x , state_type &dxdt , const double /* t */ )
{
dxdt[0] = x[1];
dxdt[1] = -x[0] - gam*x[1];
}
//]
//[ rhs_class
/* The rhs of x' = f(x) defined as a class */
class harm_osc {
double m_gam;
public:
harm_osc( double gam ) : m_gam(gam) { }
void operator() ( const state_type &x , state_type &dxdt , const double /* t */ )
{
dxdt[0] = x[1];
dxdt[1] = -x[0] - m_gam*x[1];
}
};
//]
//[ integrate_observer
struct push_back_state_and_time
{
std::vector< state_type >& m_states;
std::vector< double >& m_times;
push_back_state_and_time( std::vector< state_type > &states , std::vector< double > × )
: m_states( states ) , m_times( times ) { }
void operator()( const state_type &x , double t )
{
m_states.push_back( x );
m_times.push_back( t );
}
};
//]
struct write_state
{
void operator()( const state_type &x ) const
{
std::cout << x[0] << "\t" << x[1] << "\n";
}
};
int main(int /* argc */ , char** /* argv */ )
{
using namespace std;
using namespace boost::numeric::odeint;
//[ state_initialization
state_type x(2);
x[0] = 1.0; // start at x=1.0, p=0.0
x[1] = 0.0;
//]
//[ integration
size_t steps = integrate( harmonic_oscillator ,
x , 0.0 , 10.0 , 0.1 );
//]
//[ integration_class
harm_osc ho(0.15);
steps = integrate( ho ,
x , 0.0 , 10.0 , 0.1 );
//]
//[ integrate_observ
vector<state_type> x_vec;
vector<double> times;
steps = integrate( harmonic_oscillator ,
x , 0.0 , 10.0 , 0.1 ,
push_back_state_and_time( x_vec , times ) );
/* output */
for( size_t i=0; i<=steps; i++ )
{
cout << times[i] << '\t' << x_vec[i][0] << '\t' << x_vec[i][1] << '\n';
}
//]
//[ define_const_stepper
runge_kutta4< state_type > stepper;
integrate_const( stepper , harmonic_oscillator , x , 0.0 , 10.0 , 0.01 );
//]
//[ integrate_const_loop
const double dt = 0.01;
for( double t=0.0 ; t<10.0 ; t+= dt )
stepper.do_step( harmonic_oscillator , x , t , dt );
//]
//[ define_adapt_stepper
typedef runge_kutta_cash_karp54< state_type > error_stepper_type;
//]
//[ integrate_adapt
typedef controlled_runge_kutta< error_stepper_type > controlled_stepper_type;
controlled_stepper_type controlled_stepper;
integrate_adaptive( controlled_stepper , harmonic_oscillator , x , 0.0 , 10.0 , 0.01 );
//]
{
//[integrate_adapt_full
double abs_err = 1.0e-10 , rel_err = 1.0e-6 , a_x = 1.0 , a_dxdt = 1.0;
controlled_stepper_type controlled_stepper(
default_error_checker< double , range_algebra , default_operations >( abs_err , rel_err , a_x , a_dxdt ) );
integrate_adaptive( controlled_stepper , harmonic_oscillator , x , 0.0 , 10.0 , 0.01 );
//]
}
//[integrate_adapt_make_controlled
integrate_adaptive( make_controlled< error_stepper_type >( 1.0e-10 , 1.0e-6 ) ,
harmonic_oscillator , x , 0.0 , 10.0 , 0.01 );
//]
//[integrate_adapt_make_controlled_alternative
integrate_adaptive( make_controlled( 1.0e-10 , 1.0e-6 , error_stepper_type() ) ,
harmonic_oscillator , x , 0.0 , 10.0 , 0.01 );
//]
#ifdef BOOST_NUMERIC_ODEINT_CXX11
//[ define_const_stepper_cpp11
{
runge_kutta4< state_type > stepper;
integrate_const( stepper , []( const state_type &x , state_type &dxdt , double t ) {
dxdt[0] = x[1]; dxdt[1] = -x[0] - gam*x[1]; }
, x , 0.0 , 10.0 , 0.01 );
}
//]
//[ harm_iterator_const_step]
std::for_each( make_const_step_time_iterator_begin( stepper , harmonic_oscillator, x , 0.0 , 0.1 , 10.0 ) ,
make_const_step_time_iterator_end( stepper , harmonic_oscillator, x ) ,
[]( std::pair< const state_type & , const double & > x ) {
cout << x.second << " " << x.first[0] << " " << x.first[1] << "\n"; } );
//]
#endif
}