1<html> 2<head> 3<meta http-equiv="Content-Type" content="text/html; charset=UTF-8"> 4<title>Harmonic oscillator</title> 5<link rel="stylesheet" href="../../../../../../../doc/src/boostbook.css" type="text/css"> 6<meta name="generator" content="DocBook XSL Stylesheets V1.79.1"> 7<link rel="home" href="../../index.html" title="Chapter 1. Boost.Numeric.Odeint"> 8<link rel="up" href="../tutorial.html" title="Tutorial"> 9<link rel="prev" href="../tutorial.html" title="Tutorial"> 10<link rel="next" href="solar_system.html" title="Solar system"> 11</head> 12<body bgcolor="white" text="black" link="#0000FF" vlink="#840084" alink="#0000FF"> 13<table cellpadding="2" width="100%"><tr> 14<td valign="top"><img alt="Boost C++ Libraries" width="277" height="86" src="../../logo.jpg"></td> 15<td align="center"><a href="../../../../../../../index.html">Home</a></td> 16<td align="center"><a href="../../../../../../../libs/libraries.htm">Libraries</a></td> 17<td align="center"><a href="http://www.boost.org/users/people.html">People</a></td> 18<td align="center"><a href="http://www.boost.org/users/faq.html">FAQ</a></td> 19<td align="center"><a href="../../../../../../../more/index.htm">More</a></td> 20</tr></table> 21<hr> 22<div class="spirit-nav"> 23<a accesskey="p" href="../tutorial.html"><img src="../../../../../../../doc/src/images/prev.png" alt="Prev"></a><a accesskey="u" href="../tutorial.html"><img src="../../../../../../../doc/src/images/up.png" alt="Up"></a><a accesskey="h" href="../../index.html"><img src="../../../../../../../doc/src/images/home.png" alt="Home"></a><a accesskey="n" href="solar_system.html"><img src="../../../../../../../doc/src/images/next.png" alt="Next"></a> 24</div> 25<div class="section"> 26<div class="titlepage"><div><div><h3 class="title"> 27<a name="boost_numeric_odeint.tutorial.harmonic_oscillator"></a><a class="link" href="harmonic_oscillator.html" title="Harmonic oscillator">Harmonic 28 oscillator</a> 29</h3></div></div></div> 30<div class="toc"><dl class="toc"> 31<dt><span class="section"><a href="harmonic_oscillator.html#boost_numeric_odeint.tutorial.harmonic_oscillator.define_the_ode">Define 32 the ODE</a></span></dt> 33<dt><span class="section"><a href="harmonic_oscillator.html#boost_numeric_odeint.tutorial.harmonic_oscillator.stepper_types">Stepper 34 Types</a></span></dt> 35<dt><span class="section"><a href="harmonic_oscillator.html#boost_numeric_odeint.tutorial.harmonic_oscillator.integration_with_constant_step_size">Integration 36 with Constant Step Size</a></span></dt> 37<dt><span class="section"><a href="harmonic_oscillator.html#boost_numeric_odeint.tutorial.harmonic_oscillator.integration_with_adaptive_step_size">Integration 38 with Adaptive Step Size</a></span></dt> 39<dt><span class="section"><a href="harmonic_oscillator.html#boost_numeric_odeint.tutorial.harmonic_oscillator.using_iterators">Using 40 iterators</a></span></dt> 41</dl></div> 42<div class="section"> 43<div class="titlepage"><div><div><h4 class="title"> 44<a name="boost_numeric_odeint.tutorial.harmonic_oscillator.define_the_ode"></a><a class="link" href="harmonic_oscillator.html#boost_numeric_odeint.tutorial.harmonic_oscillator.define_the_ode" title="Define the ODE">Define 45 the ODE</a> 46</h4></div></div></div> 47<p> 48 First of all, you have to specify the data type that represents a state 49 <span class="emphasis"><em>x</em></span> of your system. Mathematically, this usually is 50 an n-dimensional vector with real numbers or complex numbers as scalar 51 objects. For odeint the most natural way is to use <code class="computeroutput"><span class="identifier">vector</span><span class="special"><</span> <span class="keyword">double</span> <span class="special">></span></code> or <code class="computeroutput"><span class="identifier">vector</span><span class="special"><</span> <span class="identifier">complex</span><span class="special"><</span> <span class="keyword">double</span> <span class="special">></span> <span class="special">></span></code> 52 to represent the system state. However, odeint can deal with other container 53 types as well, e.g. <code class="computeroutput"><span class="identifier">boost</span><span class="special">::</span><span class="identifier">array</span><span class="special"><</span> <span class="keyword">double</span> <span class="special">,</span> <span class="identifier">N</span> <span class="special">></span></code>, as long as it fulfills some requirements 54 defined below. 55 </p> 56<p> 57 To integrate a differential equation numerically, one also has to define 58 the rhs of the equation <span class="emphasis"><em>x' = f(x)</em></span>. In odeint you supply 59 this function in terms of an object that implements the ()-operator with 60 a certain parameter structure. Hence, the straightforward way would be 61 to just define a function, e.g: 62 </p> 63<p> 64</p> 65<pre class="programlisting"><span class="comment">/* The type of container used to hold the state vector */</span> 66<span class="keyword">typedef</span> <span class="identifier">std</span><span class="special">::</span><span class="identifier">vector</span><span class="special"><</span> <span class="keyword">double</span> <span class="special">></span> <span class="identifier">state_type</span><span class="special">;</span> 67 68<span class="keyword">const</span> <span class="keyword">double</span> <span class="identifier">gam</span> <span class="special">=</span> <span class="number">0.15</span><span class="special">;</span> 69 70<span class="comment">/* The rhs of x' = f(x) */</span> 71<span class="keyword">void</span> <span class="identifier">harmonic_oscillator</span><span class="special">(</span> <span class="keyword">const</span> <span class="identifier">state_type</span> <span class="special">&</span><span class="identifier">x</span> <span class="special">,</span> <span class="identifier">state_type</span> <span class="special">&</span><span class="identifier">dxdt</span> <span class="special">,</span> <span class="keyword">const</span> <span class="keyword">double</span> <span class="comment">/* t */</span> <span class="special">)</span> 72<span class="special">{</span> 73 <span class="identifier">dxdt</span><span class="special">[</span><span class="number">0</span><span class="special">]</span> <span class="special">=</span> <span class="identifier">x</span><span class="special">[</span><span class="number">1</span><span class="special">];</span> 74 <span class="identifier">dxdt</span><span class="special">[</span><span class="number">1</span><span class="special">]</span> <span class="special">=</span> <span class="special">-</span><span class="identifier">x</span><span class="special">[</span><span class="number">0</span><span class="special">]</span> <span class="special">-</span> <span class="identifier">gam</span><span class="special">*</span><span class="identifier">x</span><span class="special">[</span><span class="number">1</span><span class="special">];</span> 75<span class="special">}</span> 76</pre> 77<p> 78 </p> 79<p> 80 The parameters of the function must follow the example above where <code class="computeroutput"><span class="identifier">x</span></code> is the current state, here a two-component 81 vector containing position <span class="emphasis"><em>q</em></span> and momentum <span class="emphasis"><em>p</em></span> 82 of the oscillator, <code class="computeroutput"><span class="identifier">dxdt</span></code> 83 is the derivative <span class="emphasis"><em>x'</em></span> and should be filled by the function 84 with <span class="emphasis"><em>f(x)</em></span>, and <code class="computeroutput"><span class="identifier">t</span></code> 85 is the current time. Note that in this example <span class="emphasis"><em>t</em></span> is 86 not required to calculate <span class="emphasis"><em>f</em></span>, however odeint expects 87 the function signature to have exactly three parameters (there are exception, 88 discussed later). 89 </p> 90<p> 91 A more sophisticated approach is to implement the system as a class where 92 the rhs function is defined as the ()-operator of the class with the same 93 parameter structure as above: 94 </p> 95<p> 96</p> 97<pre class="programlisting"><span class="comment">/* The rhs of x' = f(x) defined as a class */</span> 98<span class="keyword">class</span> <span class="identifier">harm_osc</span> <span class="special">{</span> 99 100 <span class="keyword">double</span> <span class="identifier">m_gam</span><span class="special">;</span> 101 102<span class="keyword">public</span><span class="special">:</span> 103 <span class="identifier">harm_osc</span><span class="special">(</span> <span class="keyword">double</span> <span class="identifier">gam</span> <span class="special">)</span> <span class="special">:</span> <span class="identifier">m_gam</span><span class="special">(</span><span class="identifier">gam</span><span class="special">)</span> <span class="special">{</span> <span class="special">}</span> 104 105 <span class="keyword">void</span> <span class="keyword">operator</span><span class="special">()</span> <span class="special">(</span> <span class="keyword">const</span> <span class="identifier">state_type</span> <span class="special">&</span><span class="identifier">x</span> <span class="special">,</span> <span class="identifier">state_type</span> <span class="special">&</span><span class="identifier">dxdt</span> <span class="special">,</span> <span class="keyword">const</span> <span class="keyword">double</span> <span class="comment">/* t */</span> <span class="special">)</span> 106 <span class="special">{</span> 107 <span class="identifier">dxdt</span><span class="special">[</span><span class="number">0</span><span class="special">]</span> <span class="special">=</span> <span class="identifier">x</span><span class="special">[</span><span class="number">1</span><span class="special">];</span> 108 <span class="identifier">dxdt</span><span class="special">[</span><span class="number">1</span><span class="special">]</span> <span class="special">=</span> <span class="special">-</span><span class="identifier">x</span><span class="special">[</span><span class="number">0</span><span class="special">]</span> <span class="special">-</span> <span class="identifier">m_gam</span><span class="special">*</span><span class="identifier">x</span><span class="special">[</span><span class="number">1</span><span class="special">];</span> 109 <span class="special">}</span> 110<span class="special">};</span> 111</pre> 112<p> 113 </p> 114<p> 115 odeint can deal with instances of such classes instead of pure functions 116 which allows for cleaner code. 117 </p> 118</div> 119<div class="section"> 120<div class="titlepage"><div><div><h4 class="title"> 121<a name="boost_numeric_odeint.tutorial.harmonic_oscillator.stepper_types"></a><a class="link" href="harmonic_oscillator.html#boost_numeric_odeint.tutorial.harmonic_oscillator.stepper_types" title="Stepper Types">Stepper 122 Types</a> 123</h4></div></div></div> 124<p> 125 Numerical integration works iteratively, that means you start at a state 126 <span class="emphasis"><em>x(t)</em></span> and perform a time-step of length <span class="emphasis"><em>dt</em></span> 127 to obtain the approximate state <span class="emphasis"><em>x(t+dt)</em></span>. There exist 128 many different methods to perform such a time-step each of which has a 129 certain order <span class="emphasis"><em>q</em></span>. If the order of a method is <span class="emphasis"><em>q</em></span> 130 than it is accurate up to term <span class="emphasis"><em>~dt<sup>q</sup></em></span> that means the 131 error in <span class="emphasis"><em>x</em></span> made by such a step is <span class="emphasis"><em>~dt<sup>q+1</sup></em></span>. 132 odeint provides several steppers of different orders, see <a class="link" href="../odeint_in_detail/steppers.html#boost_numeric_odeint.odeint_in_detail.steppers.stepper_overview" title="Stepper overview">Stepper 133 overview</a>. 134 </p> 135<p> 136 Some of steppers in the table above are special: Some need the Jacobian 137 of the ODE, others are constructed for special ODE-systems like Hamiltonian 138 systems. We will show typical examples and use-cases in this tutorial and 139 which kind of steppers should be applied. 140 </p> 141</div> 142<div class="section"> 143<div class="titlepage"><div><div><h4 class="title"> 144<a name="boost_numeric_odeint.tutorial.harmonic_oscillator.integration_with_constant_step_size"></a><a class="link" href="harmonic_oscillator.html#boost_numeric_odeint.tutorial.harmonic_oscillator.integration_with_constant_step_size" title="Integration with Constant Step Size">Integration 145 with Constant Step Size</a> 146</h4></div></div></div> 147<p> 148 The basic stepper just performs one time-step and doesn't give you any 149 information about the error that was made (except that you know it is of 150 order <span class="emphasis"><em>q+1</em></span>). Such steppers are used with constant step 151 size that should be chosen small enough to have reasonable small errors. 152 However, you should apply some sort of validity check of your results (like 153 observing conserved quantities) because you have no other control of the 154 error. The following example defines a basic stepper based on the classical 155 Runge-Kutta scheme of 4th order. The declaration of the stepper requires 156 the state type as template parameter. The integration can now be done by 157 using the <code class="computeroutput"><span class="identifier">integrate_const</span><span class="special">(</span> <span class="identifier">Stepper</span><span class="special">,</span> <span class="identifier">System</span><span class="special">,</span> <span class="identifier">state</span><span class="special">,</span> <span class="identifier">start_time</span><span class="special">,</span> <span class="identifier">end_time</span><span class="special">,</span> <span class="identifier">step_size</span> 158 <span class="special">)</span></code> function from odeint: 159 </p> 160<p> 161</p> 162<pre class="programlisting"><span class="identifier">runge_kutta4</span><span class="special"><</span> <span class="identifier">state_type</span> <span class="special">></span> <span class="identifier">stepper</span><span class="special">;</span> 163<span class="identifier">integrate_const</span><span class="special">(</span> <span class="identifier">stepper</span> <span class="special">,</span> <span class="identifier">harmonic_oscillator</span> <span class="special">,</span> <span class="identifier">x</span> <span class="special">,</span> <span class="number">0.0</span> <span class="special">,</span> <span class="number">10.0</span> <span class="special">,</span> <span class="number">0.01</span> <span class="special">);</span> 164</pre> 165<p> 166 </p> 167<p> 168 This call integrates the system defined by <code class="computeroutput"><span class="identifier">harmonic_oscillator</span></code> 169 using the RK4 method from <span class="emphasis"><em>t=0</em></span> to <span class="emphasis"><em>10</em></span> 170 with a step-size <span class="emphasis"><em>dt=0.01</em></span> and the initial condition 171 given in <code class="computeroutput"><span class="identifier">x</span></code>. The result, 172 <span class="emphasis"><em>x(t=10)</em></span> is stored in <code class="computeroutput"><span class="identifier">x</span></code> 173 (in-place). Each stepper defines a <code class="computeroutput"><span class="identifier">do_step</span></code> 174 method which can also be used directly. So, you write down the above example 175 as 176 </p> 177<p> 178</p> 179<pre class="programlisting"><span class="keyword">const</span> <span class="keyword">double</span> <span class="identifier">dt</span> <span class="special">=</span> <span class="number">0.01</span><span class="special">;</span> 180<span class="keyword">for</span><span class="special">(</span> <span class="keyword">double</span> <span class="identifier">t</span><span class="special">=</span><span class="number">0.0</span> <span class="special">;</span> <span class="identifier">t</span><span class="special"><</span><span class="number">10.0</span> <span class="special">;</span> <span class="identifier">t</span><span class="special">+=</span> <span class="identifier">dt</span> <span class="special">)</span> 181 <span class="identifier">stepper</span><span class="special">.</span><span class="identifier">do_step</span><span class="special">(</span> <span class="identifier">harmonic_oscillator</span> <span class="special">,</span> <span class="identifier">x</span> <span class="special">,</span> <span class="identifier">t</span> <span class="special">,</span> <span class="identifier">dt</span> <span class="special">);</span> 182</pre> 183<p> 184 </p> 185<div class="tip"><table border="0" summary="Tip"> 186<tr> 187<td rowspan="2" align="center" valign="top" width="25"><img alt="[Tip]" src="../../../../../../../doc/src/images/tip.png"></td> 188<th align="left">Tip</th> 189</tr> 190<tr><td align="left" valign="top"> 191<p> 192 If you have a C++11 enabled compiler you can easily use lambdas to create 193 the system function : 194 </p> 195<p> 196</p> 197<pre class="programlisting"><span class="special">{</span> 198<span class="identifier">runge_kutta4</span><span class="special"><</span> <span class="identifier">state_type</span> <span class="special">></span> <span class="identifier">stepper</span><span class="special">;</span> 199<span class="identifier">integrate_const</span><span class="special">(</span> <span class="identifier">stepper</span> <span class="special">,</span> <span class="special">[](</span> <span class="keyword">const</span> <span class="identifier">state_type</span> <span class="special">&</span><span class="identifier">x</span> <span class="special">,</span> <span class="identifier">state_type</span> <span class="special">&</span><span class="identifier">dxdt</span> <span class="special">,</span> <span class="keyword">double</span> <span class="identifier">t</span> <span class="special">)</span> <span class="special">{</span> 200 <span class="identifier">dxdt</span><span class="special">[</span><span class="number">0</span><span class="special">]</span> <span class="special">=</span> <span class="identifier">x</span><span class="special">[</span><span class="number">1</span><span class="special">];</span> <span class="identifier">dxdt</span><span class="special">[</span><span class="number">1</span><span class="special">]</span> <span class="special">=</span> <span class="special">-</span><span class="identifier">x</span><span class="special">[</span><span class="number">0</span><span class="special">]</span> <span class="special">-</span> <span class="identifier">gam</span><span class="special">*</span><span class="identifier">x</span><span class="special">[</span><span class="number">1</span><span class="special">];</span> <span class="special">}</span> 201 <span class="special">,</span> <span class="identifier">x</span> <span class="special">,</span> <span class="number">0.0</span> <span class="special">,</span> <span class="number">10.0</span> <span class="special">,</span> <span class="number">0.01</span> <span class="special">);</span> 202<span class="special">}</span> 203</pre> 204<p> 205 </p> 206</td></tr> 207</table></div> 208</div> 209<div class="section"> 210<div class="titlepage"><div><div><h4 class="title"> 211<a name="boost_numeric_odeint.tutorial.harmonic_oscillator.integration_with_adaptive_step_size"></a><a class="link" href="harmonic_oscillator.html#boost_numeric_odeint.tutorial.harmonic_oscillator.integration_with_adaptive_step_size" title="Integration with Adaptive Step Size">Integration 212 with Adaptive Step Size</a> 213</h4></div></div></div> 214<p> 215 To improve the numerical results and additionally minimize the computational 216 effort, the application of a step size control is advisable. Step size 217 control is realized via stepper algorithms that additionally provide an 218 error estimation of the applied step. odeint provides a number of such 219 <span class="bold"><strong>ErrorSteppers</strong></span> and we will show their usage 220 on the example of <code class="computeroutput"><span class="identifier">explicit_error_rk54_ck</span></code> 221 - a 5th order Runge-Kutta method with 4th order error estimation and coefficients 222 introduced by Cash and Karp. 223 </p> 224<p> 225</p> 226<pre class="programlisting"><span class="keyword">typedef</span> <span class="identifier">runge_kutta_cash_karp54</span><span class="special"><</span> <span class="identifier">state_type</span> <span class="special">></span> <span class="identifier">error_stepper_type</span><span class="special">;</span> 227</pre> 228<p> 229 </p> 230<p> 231 Given the error stepper, one still needs an instance that checks the error 232 and adjusts the step size accordingly. In odeint, this is done by <span class="bold"><strong>ControlledSteppers</strong></span>. For the <code class="computeroutput"><span class="identifier">runge_kutta_cash_karp54</span></code> 233 stepper a <code class="computeroutput"><span class="identifier">controlled_runge_kutta</span></code> 234 stepper exists which can be used via 235 </p> 236<p> 237</p> 238<pre class="programlisting"><span class="keyword">typedef</span> <span class="identifier">controlled_runge_kutta</span><span class="special"><</span> <span class="identifier">error_stepper_type</span> <span class="special">></span> <span class="identifier">controlled_stepper_type</span><span class="special">;</span> 239<span class="identifier">controlled_stepper_type</span> <span class="identifier">controlled_stepper</span><span class="special">;</span> 240<span class="identifier">integrate_adaptive</span><span class="special">(</span> <span class="identifier">controlled_stepper</span> <span class="special">,</span> <span class="identifier">harmonic_oscillator</span> <span class="special">,</span> <span class="identifier">x</span> <span class="special">,</span> <span class="number">0.0</span> <span class="special">,</span> <span class="number">10.0</span> <span class="special">,</span> <span class="number">0.01</span> <span class="special">);</span> 241</pre> 242<p> 243 </p> 244<p> 245 As above, this integrates the system defined by <code class="computeroutput"><span class="identifier">harmonic_oscillator</span></code>, 246 but now using an adaptive step size method based on the Runge-Kutta Cash-Karp 247 54 scheme from <span class="emphasis"><em>t=0</em></span> to <span class="emphasis"><em>10</em></span> with 248 an initial step size of <span class="emphasis"><em>dt=0.01</em></span> (will be adjusted) 249 and the initial condition given in x. The result, <span class="emphasis"><em>x(t=10)</em></span>, 250 will also be stored in x (in-place). 251 </p> 252<p> 253 In the above example an error stepper is nested in a controlled stepper. 254 This is a nice technique; however one drawback is that one always needs 255 to define both steppers. One could also write the instantiation of the 256 controlled stepper into the call of the integrate function but a complete 257 knowledge of the underlying stepper types is still necessary. Another point 258 is, that the error tolerances for the step size control are not easily 259 included into the controlled stepper. Both issues can be solved by using 260 <code class="computeroutput"><span class="identifier">make_controlled</span></code>: 261 </p> 262<p> 263</p> 264<pre class="programlisting"><span class="identifier">integrate_adaptive</span><span class="special">(</span> <span class="identifier">make_controlled</span><span class="special"><</span> <span class="identifier">error_stepper_type</span> <span class="special">>(</span> <span class="number">1.0e-10</span> <span class="special">,</span> <span class="number">1.0e-6</span> <span class="special">)</span> <span class="special">,</span> 265 <span class="identifier">harmonic_oscillator</span> <span class="special">,</span> <span class="identifier">x</span> <span class="special">,</span> <span class="number">0.0</span> <span class="special">,</span> <span class="number">10.0</span> <span class="special">,</span> <span class="number">0.01</span> <span class="special">);</span> 266</pre> 267<p> 268 </p> 269<p> 270 <code class="computeroutput"><span class="identifier">make_controlled</span></code> can be 271 used with many of the steppers of odeint. The first parameter is the absolute 272 error tolerance <span class="emphasis"><em>eps_abs</em></span> and the second is the relative 273 error tolerance <span class="emphasis"><em>eps_rel</em></span> which is used during the integration. 274 The template parameter determines from which error stepper a controlled 275 stepper should be instantiated. An alternative syntax of <code class="computeroutput"><span class="identifier">make_controlled</span></code> is 276 </p> 277<p> 278</p> 279<pre class="programlisting"><span class="identifier">integrate_adaptive</span><span class="special">(</span> <span class="identifier">make_controlled</span><span class="special">(</span> <span class="number">1.0e-10</span> <span class="special">,</span> <span class="number">1.0e-6</span> <span class="special">,</span> <span class="identifier">error_stepper_type</span><span class="special">()</span> <span class="special">)</span> <span class="special">,</span> 280 <span class="identifier">harmonic_oscillator</span> <span class="special">,</span> <span class="identifier">x</span> <span class="special">,</span> <span class="number">0.0</span> <span class="special">,</span> <span class="number">10.0</span> <span class="special">,</span> <span class="number">0.01</span> <span class="special">);</span> 281</pre> 282<p> 283 </p> 284<p> 285 For the Runge-Kutta controller the error made during one step is compared 286 with <span class="emphasis"><em>eps_abs + eps_rel * ( a<sub>x</sub> * |x| + a<sub>dxdt</sub> * dt * |dxdt| )</em></span>. 287 If the error is smaller than this value the current step is accepted, otherwise 288 it is rejected and the step size is decreased. Note, that the step size 289 is also increased if the error gets too small compared to the rhs of the 290 above relation. The full instantiation of the <code class="computeroutput"><span class="identifier">controlled_runge_kutta</span></code> 291 with all parameters is therefore 292 </p> 293<p> 294</p> 295<pre class="programlisting"><span class="keyword">double</span> <span class="identifier">abs_err</span> <span class="special">=</span> <span class="number">1.0e-10</span> <span class="special">,</span> <span class="identifier">rel_err</span> <span class="special">=</span> <span class="number">1.0e-6</span> <span class="special">,</span> <span class="identifier">a_x</span> <span class="special">=</span> <span class="number">1.0</span> <span class="special">,</span> <span class="identifier">a_dxdt</span> <span class="special">=</span> <span class="number">1.0</span><span class="special">;</span> 296<span class="identifier">controlled_stepper_type</span> <span class="identifier">controlled_stepper</span><span class="special">(</span> 297 <span class="identifier">default_error_checker</span><span class="special"><</span> <span class="keyword">double</span> <span class="special">,</span> <span class="identifier">range_algebra</span> <span class="special">,</span> <span class="identifier">default_operations</span> <span class="special">>(</span> <span class="identifier">abs_err</span> <span class="special">,</span> <span class="identifier">rel_err</span> <span class="special">,</span> <span class="identifier">a_x</span> <span class="special">,</span> <span class="identifier">a_dxdt</span> <span class="special">)</span> <span class="special">);</span> 298<span class="identifier">integrate_adaptive</span><span class="special">(</span> <span class="identifier">controlled_stepper</span> <span class="special">,</span> <span class="identifier">harmonic_oscillator</span> <span class="special">,</span> <span class="identifier">x</span> <span class="special">,</span> <span class="number">0.0</span> <span class="special">,</span> <span class="number">10.0</span> <span class="special">,</span> <span class="number">0.01</span> <span class="special">);</span> 299</pre> 300<p> 301 </p> 302<p> 303 When using <code class="computeroutput"><span class="identifier">make_controlled</span></code> 304 the parameter <span class="emphasis"><em>a<sub>x</sub></em></span> and <span class="emphasis"><em>a<sub>dxdt</sub></em></span> are 305 used with their standard values of 1. 306 </p> 307<p> 308 In the tables below, one can find all steppers which are working with 309 <code class="computeroutput"><span class="identifier">make_controlled</span></code> and <code class="computeroutput"><span class="identifier">make_dense_output</span></code> which is the analog 310 for the dense output steppers. 311 </p> 312<div class="table"> 313<a name="boost_numeric_odeint.tutorial.harmonic_oscillator.integration_with_adaptive_step_size.generation_functions_make_controlled__abs_error___rel_error___stepper__"></a><p class="title"><b>Table 1.2. Generation functions make_controlled( abs_error , rel_error , stepper 314 )</b></p> 315<div class="table-contents"><table class="table" summary="Generation functions make_controlled( abs_error , rel_error , stepper 316 )"> 317<colgroup> 318<col> 319<col> 320<col> 321</colgroup> 322<thead><tr> 323<th> 324 <p> 325 Stepper 326 </p> 327 </th> 328<th> 329 <p> 330 Result of make_controlled 331 </p> 332 </th> 333<th> 334 <p> 335 Remarks 336 </p> 337 </th> 338</tr></thead> 339<tbody> 340<tr> 341<td> 342 <p> 343 <code class="computeroutput"><span class="identifier">runge_kutta_cash_karp54</span></code> 344 </p> 345 </td> 346<td> 347 <p> 348 <code class="computeroutput"><span class="identifier">controlled_runge_kutta</span><span class="special"><</span> <span class="identifier">runge_kutta_cash_karp54</span> 349 <span class="special">,</span> <span class="identifier">default_error_checker</span><span class="special"><...></span> <span class="special">></span></code> 350 </p> 351 </td> 352<td> 353 <p> 354 <span class="emphasis"><em>a<sub>x</sub>=1</em></span>, <span class="emphasis"><em>a<sub>dxdt</sub>=1</em></span> 355 </p> 356 </td> 357</tr> 358<tr> 359<td> 360 <p> 361 <code class="computeroutput"><span class="identifier">runge_kutta_fehlberg78</span></code> 362 </p> 363 </td> 364<td> 365 <p> 366 <code class="computeroutput"><span class="identifier">controlled_runge_kutta</span><span class="special"><</span> <span class="identifier">runge_kutta_fehlberg78</span> 367 <span class="special">,</span> <span class="identifier">default_error_checker</span><span class="special"><...></span> <span class="special">></span></code> 368 </p> 369 </td> 370<td> 371 <p> 372 <span class="emphasis"><em>a<sub>x</sub>=1</em></span>, <span class="emphasis"><em>a<sub>dxdt</sub>=1</em></span> 373 </p> 374 </td> 375</tr> 376<tr> 377<td> 378 <p> 379 <code class="computeroutput"><span class="identifier">runge_kutta_dopri5</span></code> 380 </p> 381 </td> 382<td> 383 <p> 384 <code class="computeroutput"><span class="identifier">controlled_runge_kutta</span><span class="special"><</span> <span class="identifier">runge_kutta_dopri5</span> 385 <span class="special">,</span> <span class="identifier">default_error_checker</span><span class="special"><...></span> <span class="special">></span></code> 386 </p> 387 </td> 388<td> 389 <p> 390 <span class="emphasis"><em>a <sub>x</sub>=1</em></span>, <span class="emphasis"><em>a<sub>dxdt</sub>=1</em></span> 391 </p> 392 </td> 393</tr> 394<tr> 395<td> 396 <p> 397 <code class="computeroutput"><span class="identifier">rosenbrock4</span></code> 398 </p> 399 </td> 400<td> 401 <p> 402 <code class="computeroutput"><span class="identifier">rosenbrock4_controlled</span><span class="special"><</span> <span class="identifier">rosenbrock4</span> 403 <span class="special">></span></code> 404 </p> 405 </td> 406<td> 407 <p> 408 - 409 </p> 410 </td> 411</tr> 412</tbody> 413</table></div> 414</div> 415<br class="table-break"><div class="table"> 416<a name="boost_numeric_odeint.tutorial.harmonic_oscillator.integration_with_adaptive_step_size.generation_functions_make_dense_output__abs_error___rel_error___stepper__"></a><p class="title"><b>Table 1.3. Generation functions make_dense_output( abs_error , rel_error , 417 stepper )</b></p> 418<div class="table-contents"><table class="table" summary="Generation functions make_dense_output( abs_error , rel_error , 419 stepper )"> 420<colgroup> 421<col> 422<col> 423<col> 424</colgroup> 425<thead><tr> 426<th> 427 <p> 428 Stepper 429 </p> 430 </th> 431<th> 432 <p> 433 Result of make_dense_output 434 </p> 435 </th> 436<th> 437 <p> 438 Remarks 439 </p> 440 </th> 441</tr></thead> 442<tbody> 443<tr> 444<td> 445 <p> 446 <code class="computeroutput"><span class="identifier">runge_kutta_dopri5</span></code> 447 </p> 448 </td> 449<td> 450 <p> 451 <code class="computeroutput"><span class="identifier">dense_output_runge_kutta</span><span class="special"><</span> <span class="identifier">controlled_runge_kutta</span><span class="special"><</span> <span class="identifier">runge_kutta_dopri5</span> 452 <span class="special">,</span> <span class="identifier">default_error_checker</span><span class="special"><...></span> <span class="special">></span> 453 <span class="special">></span></code> 454 </p> 455 </td> 456<td> 457 <p> 458 <span class="emphasis"><em>a <sub>x</sub>=1</em></span>, <span class="emphasis"><em>a<sub>dxdt</sub>=1</em></span> 459 </p> 460 </td> 461</tr> 462<tr> 463<td> 464 <p> 465 <code class="computeroutput"><span class="identifier">rosenbrock4</span></code> 466 </p> 467 </td> 468<td> 469 <p> 470 <code class="computeroutput"><span class="identifier">rosenbrock4_dense_output</span><span class="special"><</span> <span class="identifier">rosenbrock4_controller</span><span class="special"><</span> <span class="identifier">rosenbrock4</span> 471 <span class="special">></span> <span class="special">></span></code> 472 </p> 473 </td> 474<td> 475 <p> 476 - 477 </p> 478 </td> 479</tr> 480</tbody> 481</table></div> 482</div> 483<br class="table-break"><p> 484 When using <code class="computeroutput"><span class="identifier">make_controlled</span></code> 485 or <code class="computeroutput"><span class="identifier">make_dense_output</span></code> one 486 should be aware which exact type is used and how the step size control 487 works. 488 </p> 489</div> 490<div class="section"> 491<div class="titlepage"><div><div><h4 class="title"> 492<a name="boost_numeric_odeint.tutorial.harmonic_oscillator.using_iterators"></a><a class="link" href="harmonic_oscillator.html#boost_numeric_odeint.tutorial.harmonic_oscillator.using_iterators" title="Using iterators">Using 493 iterators</a> 494</h4></div></div></div> 495<p> 496 odeint supports iterators for solving ODEs. That is, you instantiate a 497 pair of iterators and instead of using the integrate routines with an appropriate 498 observer you put the iterators in one of the algorithm from the C++ standard 499 library or from Boost.Range. An example is 500 </p> 501<p> 502</p> 503<pre class="programlisting"><span class="identifier">std</span><span class="special">::</span><span class="identifier">for_each</span><span class="special">(</span> <span class="identifier">make_const_step_time_iterator_begin</span><span class="special">(</span> <span class="identifier">stepper</span> <span class="special">,</span> <span class="identifier">harmonic_oscillator</span><span class="special">,</span> <span class="identifier">x</span> <span class="special">,</span> <span class="number">0.0</span> <span class="special">,</span> <span class="number">0.1</span> <span class="special">,</span> <span class="number">10.0</span> <span class="special">)</span> <span class="special">,</span> 504 <span class="identifier">make_const_step_time_iterator_end</span><span class="special">(</span> <span class="identifier">stepper</span> <span class="special">,</span> <span class="identifier">harmonic_oscillator</span><span class="special">,</span> <span class="identifier">x</span> <span class="special">)</span> <span class="special">,</span> 505 <span class="special">[](</span> <span class="identifier">std</span><span class="special">::</span><span class="identifier">pair</span><span class="special"><</span> <span class="keyword">const</span> <span class="identifier">state_type</span> <span class="special">&</span> <span class="special">,</span> <span class="keyword">const</span> <span class="keyword">double</span> <span class="special">&</span> <span class="special">></span> <span class="identifier">x</span> <span class="special">)</span> <span class="special">{</span> 506 <span class="identifier">cout</span> <span class="special"><<</span> <span class="identifier">x</span><span class="special">.</span><span class="identifier">second</span> <span class="special"><<</span> <span class="string">" "</span> <span class="special"><<</span> <span class="identifier">x</span><span class="special">.</span><span class="identifier">first</span><span class="special">[</span><span class="number">0</span><span class="special">]</span> <span class="special"><<</span> <span class="string">" "</span> <span class="special"><<</span> <span class="identifier">x</span><span class="special">.</span><span class="identifier">first</span><span class="special">[</span><span class="number">1</span><span class="special">]</span> <span class="special"><<</span> <span class="string">"\n"</span><span class="special">;</span> <span class="special">}</span> <span class="special">);</span> 507</pre> 508<p> 509 </p> 510</div> 511<p> 512 The full source file for this example can be found here: <a href="https://github.com/headmyshoulder/odeint-v2/blob/master/examples/harmonic_oscillator.cpp" target="_top">harmonic_oscillator.cpp</a> 513 </p> 514</div> 515<table xmlns:rev="http://www.cs.rpi.edu/~gregod/boost/tools/doc/revision" width="100%"><tr> 516<td align="left"></td> 517<td align="right"><div class="copyright-footer">Copyright © 2009-2015 Karsten Ahnert and Mario Mulansky<p> 518 Distributed under the Boost Software License, Version 1.0. (See accompanying 519 file LICENSE_1_0.txt or copy at <a href="http://www.boost.org/LICENSE_1_0.txt" target="_top">http://www.boost.org/LICENSE_1_0.txt</a>) 520 </p> 521</div></td> 522</tr></table> 523<hr> 524<div class="spirit-nav"> 525<a accesskey="p" href="../tutorial.html"><img src="../../../../../../../doc/src/images/prev.png" alt="Prev"></a><a accesskey="u" href="../tutorial.html"><img src="../../../../../../../doc/src/images/up.png" alt="Up"></a><a accesskey="h" href="../../index.html"><img src="../../../../../../../doc/src/images/home.png" alt="Home"></a><a accesskey="n" href="solar_system.html"><img src="../../../../../../../doc/src/images/next.png" alt="Next"></a> 526</div> 527</body> 528</html> 529