summaryrefslogtreecommitdiff
path: root/libs/numeric/odeint/doc/html/boost_numeric_odeint/tutorial/stiff_systems.html
blob: 1f68df3a6d424981c6f87a7a394e2c99a9e69524 (plain)
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
<html>
<head>
<meta http-equiv="Content-Type" content="text/html; charset=US-ASCII">
<title>Stiff systems</title>
<link rel="stylesheet" href="../../../../../../../doc/src/boostbook.css" type="text/css">
<meta name="generator" content="DocBook XSL Stylesheets V1.78.1">
<link rel="home" href="../../index.html" title="Chapter&#160;1.&#160;Boost.Numeric.Odeint">
<link rel="up" href="../tutorial.html" title="Tutorial">
<link rel="prev" href="chaotic_systems_and_lyapunov_exponents.html" title="Chaotic systems and Lyapunov exponents">
<link rel="next" href="complex_state_types.html" title="Complex state types">
</head>
<body bgcolor="white" text="black" link="#0000FF" vlink="#840084" alink="#0000FF">
<table cellpadding="2" width="100%"><tr>
<td valign="top"><img alt="Boost C++ Libraries" width="277" height="86" src="../../logo.jpg"></td>
<td align="center"><a href="../../../../../../../index.html">Home</a></td>
<td align="center"><a href="../../../../../../../libs/libraries.htm">Libraries</a></td>
<td align="center"><a href="http://www.boost.org/users/people.html">People</a></td>
<td align="center"><a href="http://www.boost.org/users/faq.html">FAQ</a></td>
<td align="center"><a href="../../../../../../../more/index.htm">More</a></td>
</tr></table>
<hr>
<div class="spirit-nav">
<a accesskey="p" href="chaotic_systems_and_lyapunov_exponents.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="complex_state_types.html"><img src="../../../../../../../doc/src/images/next.png" alt="Next"></a>
</div>
<div class="section">
<div class="titlepage"><div><div><h3 class="title">
<a name="boost_numeric_odeint.tutorial.stiff_systems"></a><a class="link" href="stiff_systems.html" title="Stiff systems">Stiff systems</a>
</h3></div></div></div>
<p>
        An important class of ordinary differential equations are so called stiff
        system which are characterized by two or more time scales of different order.
        Examples of such systems are found in chemical systems where reaction rates
        of individual sub-reaction might differ over large ranges, for example:
      </p>
<p>
        <span class="emphasis"><em>d S<sub>&#8203;1</sub> / dt = - 101 S<sub>&#8203;2</sub> - 100 S<sub>&#8203;1</sub></em></span>
      </p>
<p>
        <span class="emphasis"><em>d S<sub>&#8203;2</sub> / dt = S<sub>&#8203;1</sub></em></span>
      </p>
<p>
        In order to efficiently solve stiff systems numerically the Jacobian
      </p>
<p>
        <span class="emphasis"><em>J = d f<sub>&#8203;i</sub> / d x<sub>&#8203;j</sub></em></span>
      </p>
<p>
        is needed. Here is the definition of the above example
      </p>
<p>
</p>
<pre class="programlisting"><span class="keyword">typedef</span> <span class="identifier">boost</span><span class="special">::</span><span class="identifier">numeric</span><span class="special">::</span><span class="identifier">ublas</span><span class="special">::</span><span class="identifier">vector</span><span class="special">&lt;</span> <span class="keyword">double</span> <span class="special">&gt;</span> <span class="identifier">vector_type</span><span class="special">;</span>
<span class="keyword">typedef</span> <span class="identifier">boost</span><span class="special">::</span><span class="identifier">numeric</span><span class="special">::</span><span class="identifier">ublas</span><span class="special">::</span><span class="identifier">matrix</span><span class="special">&lt;</span> <span class="keyword">double</span> <span class="special">&gt;</span> <span class="identifier">matrix_type</span><span class="special">;</span>

<span class="keyword">struct</span> <span class="identifier">stiff_system</span>
<span class="special">{</span>
    <span class="keyword">void</span> <span class="keyword">operator</span><span class="special">()(</span> <span class="keyword">const</span> <span class="identifier">vector_type</span> <span class="special">&amp;</span><span class="identifier">x</span> <span class="special">,</span> <span class="identifier">vector_type</span> <span class="special">&amp;</span><span class="identifier">dxdt</span> <span class="special">,</span> <span class="keyword">double</span> <span class="comment">/* t */</span> <span class="special">)</span>
    <span class="special">{</span>
        <span class="identifier">dxdt</span><span class="special">[</span> <span class="number">0</span> <span class="special">]</span> <span class="special">=</span> <span class="special">-</span><span class="number">101.0</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="number">100.0</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="identifier">x</span><span class="special">[</span> <span class="number">0</span> <span class="special">];</span>
    <span class="special">}</span>
<span class="special">};</span>

<span class="keyword">struct</span> <span class="identifier">stiff_system_jacobi</span>
<span class="special">{</span>
    <span class="keyword">void</span> <span class="keyword">operator</span><span class="special">()(</span> <span class="keyword">const</span> <span class="identifier">vector_type</span> <span class="special">&amp;</span> <span class="comment">/* x */</span> <span class="special">,</span> <span class="identifier">matrix_type</span> <span class="special">&amp;</span><span class="identifier">J</span> <span class="special">,</span> <span class="keyword">const</span> <span class="keyword">double</span> <span class="special">&amp;</span> <span class="comment">/* t */</span> <span class="special">,</span> <span class="identifier">vector_type</span> <span class="special">&amp;</span><span class="identifier">dfdt</span> <span class="special">)</span>
    <span class="special">{</span>
        <span class="identifier">J</span><span class="special">(</span> <span class="number">0</span> <span class="special">,</span> <span class="number">0</span> <span class="special">)</span> <span class="special">=</span> <span class="special">-</span><span class="number">101.0</span><span class="special">;</span>
        <span class="identifier">J</span><span class="special">(</span> <span class="number">0</span> <span class="special">,</span> <span class="number">1</span> <span class="special">)</span> <span class="special">=</span> <span class="special">-</span><span class="number">100.0</span><span class="special">;</span>
        <span class="identifier">J</span><span class="special">(</span> <span class="number">1</span> <span class="special">,</span> <span class="number">0</span> <span class="special">)</span> <span class="special">=</span> <span class="number">1.0</span><span class="special">;</span>
        <span class="identifier">J</span><span class="special">(</span> <span class="number">1</span> <span class="special">,</span> <span class="number">1</span> <span class="special">)</span> <span class="special">=</span> <span class="number">0.0</span><span class="special">;</span>
        <span class="identifier">dfdt</span><span class="special">[</span><span class="number">0</span><span class="special">]</span> <span class="special">=</span> <span class="number">0.0</span><span class="special">;</span>
        <span class="identifier">dfdt</span><span class="special">[</span><span class="number">1</span><span class="special">]</span> <span class="special">=</span> <span class="number">0.0</span><span class="special">;</span>
    <span class="special">}</span>
<span class="special">};</span>
</pre>
<p>
      </p>
<p>
        The state type has to be a <code class="computeroutput"><span class="identifier">ublas</span><span class="special">::</span><span class="identifier">vector</span></code>
        and the matrix type must by a <code class="computeroutput"><span class="identifier">ublas</span><span class="special">::</span><span class="identifier">matrix</span></code>
        since the stiff integrator only accepts these types. However, you might want
        use non-stiff integrators on this system, too - we will do so later for demonstration.
        Therefore we want to use the same function also with other state_types, realized
        by templatizing the <code class="computeroutput"><span class="keyword">operator</span><span class="special">()</span></code>:
      </p>
<p>
</p>
<pre class="programlisting"><span class="keyword">typedef</span> <span class="identifier">boost</span><span class="special">::</span><span class="identifier">numeric</span><span class="special">::</span><span class="identifier">ublas</span><span class="special">::</span><span class="identifier">vector</span><span class="special">&lt;</span> <span class="keyword">double</span> <span class="special">&gt;</span> <span class="identifier">vector_type</span><span class="special">;</span>
<span class="keyword">typedef</span> <span class="identifier">boost</span><span class="special">::</span><span class="identifier">numeric</span><span class="special">::</span><span class="identifier">ublas</span><span class="special">::</span><span class="identifier">matrix</span><span class="special">&lt;</span> <span class="keyword">double</span> <span class="special">&gt;</span> <span class="identifier">matrix_type</span><span class="special">;</span>

<span class="keyword">struct</span> <span class="identifier">stiff_system</span>
<span class="special">{</span>
    <span class="keyword">template</span><span class="special">&lt;</span> <span class="keyword">class</span> <span class="identifier">State</span> <span class="special">&gt;</span>
    <span class="keyword">void</span> <span class="keyword">operator</span><span class="special">()(</span> <span class="keyword">const</span> <span class="identifier">State</span> <span class="special">&amp;</span><span class="identifier">x</span> <span class="special">,</span> <span class="identifier">State</span> <span class="special">&amp;</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>
        <span class="special">...</span>
    <span class="special">}</span>
<span class="special">};</span>

<span class="keyword">struct</span> <span class="identifier">stiff_system_jacobi</span>
<span class="special">{</span>
    <span class="keyword">template</span><span class="special">&lt;</span> <span class="keyword">class</span> <span class="identifier">State</span> <span class="special">,</span> <span class="keyword">class</span> <span class="identifier">Matrix</span> <span class="special">&gt;</span>
    <span class="keyword">void</span> <span class="keyword">operator</span><span class="special">()(</span> <span class="keyword">const</span> <span class="identifier">State</span> <span class="special">&amp;</span><span class="identifier">x</span> <span class="special">,</span> <span class="identifier">Matrix</span> <span class="special">&amp;</span><span class="identifier">J</span> <span class="special">,</span> <span class="keyword">const</span> <span class="keyword">double</span> <span class="special">&amp;</span><span class="identifier">t</span> <span class="special">,</span> <span class="identifier">State</span> <span class="special">&amp;</span><span class="identifier">dfdt</span> <span class="special">)</span>
    <span class="special">{</span>
        <span class="special">...</span>
    <span class="special">}</span>
<span class="special">};</span>
</pre>
<p>
      </p>
<p>
        Now you can use <code class="computeroutput"><span class="identifier">stiff_system</span></code>
        in combination with <code class="computeroutput"><span class="identifier">std</span><span class="special">::</span><span class="identifier">vector</span></code> or <code class="computeroutput"><span class="identifier">boost</span><span class="special">::</span><span class="identifier">array</span></code>.
        In the example the explicit time derivative of <span class="emphasis"><em>f(x,t)</em></span>
        is introduced separately in the Jacobian. If <span class="emphasis"><em>df / dt = 0</em></span>
        simply fill <code class="computeroutput"><span class="identifier">dfdt</span></code> with zeros.
      </p>
<p>
        A well know solver for stiff systems is the Rosenbrock method. It has a step
        size control and dense output facilities and can be used like all the other
        steppers:
      </p>
<p>
</p>
<pre class="programlisting"><span class="identifier">vector_type</span> <span class="identifier">x</span><span class="special">(</span> <span class="number">2</span> <span class="special">,</span> <span class="number">1.0</span> <span class="special">);</span>

<span class="identifier">size_t</span> <span class="identifier">num_of_steps</span> <span class="special">=</span> <span class="identifier">integrate_const</span><span class="special">(</span> <span class="identifier">make_dense_output</span><span class="special">&lt;</span> <span class="identifier">rosenbrock4</span><span class="special">&lt;</span> <span class="keyword">double</span> <span class="special">&gt;</span> <span class="special">&gt;(</span> <span class="number">1.0e-6</span> <span class="special">,</span> <span class="number">1.0e-6</span> <span class="special">)</span> <span class="special">,</span>
        <span class="identifier">make_pair</span><span class="special">(</span> <span class="identifier">stiff_system</span><span class="special">()</span> <span class="special">,</span> <span class="identifier">stiff_system_jacobi</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.0</span> <span class="special">,</span> <span class="number">50.0</span> <span class="special">,</span> <span class="number">0.01</span> <span class="special">,</span>
        <span class="identifier">cout</span> <span class="special">&lt;&lt;</span> <span class="identifier">phoenix</span><span class="special">::</span><span class="identifier">arg_names</span><span class="special">::</span><span class="identifier">arg2</span> <span class="special">&lt;&lt;</span> <span class="string">" "</span> <span class="special">&lt;&lt;</span> <span class="identifier">phoenix</span><span class="special">::</span><span class="identifier">arg_names</span><span class="special">::</span><span class="identifier">arg1</span><span class="special">[</span><span class="number">0</span><span class="special">]</span> <span class="special">&lt;&lt;</span> <span class="string">"\n"</span> <span class="special">);</span>
</pre>
<p>
      </p>
<p>
        During the integration 71 steps have been done. Comparing to a classical
        Runge-Kutta solver this is a very good result. For example the Dormand-Prince
        5 method with step size control and dense output yields 1531 steps.
      </p>
<p>
</p>
<pre class="programlisting"><span class="identifier">vector_type</span> <span class="identifier">x2</span><span class="special">(</span> <span class="number">2</span> <span class="special">,</span> <span class="number">1.0</span> <span class="special">);</span>

<span class="identifier">size_t</span> <span class="identifier">num_of_steps2</span> <span class="special">=</span> <span class="identifier">integrate_const</span><span class="special">(</span> <span class="identifier">make_dense_output</span><span class="special">&lt;</span> <span class="identifier">runge_kutta_dopri5</span><span class="special">&lt;</span> <span class="identifier">vector_type</span> <span class="special">&gt;</span> <span class="special">&gt;(</span> <span class="number">1.0e-6</span> <span class="special">,</span> <span class="number">1.0e-6</span> <span class="special">)</span> <span class="special">,</span>
        <span class="identifier">stiff_system</span><span class="special">()</span> <span class="special">,</span> <span class="identifier">x2</span> <span class="special">,</span> <span class="number">0.0</span> <span class="special">,</span> <span class="number">50.0</span> <span class="special">,</span> <span class="number">0.01</span> <span class="special">,</span>
        <span class="identifier">cout</span> <span class="special">&lt;&lt;</span> <span class="identifier">phoenix</span><span class="special">::</span><span class="identifier">arg_names</span><span class="special">::</span><span class="identifier">arg2</span> <span class="special">&lt;&lt;</span> <span class="string">" "</span> <span class="special">&lt;&lt;</span> <span class="identifier">phoenix</span><span class="special">::</span><span class="identifier">arg_names</span><span class="special">::</span><span class="identifier">arg1</span><span class="special">[</span><span class="number">0</span><span class="special">]</span> <span class="special">&lt;&lt;</span> <span class="string">"\n"</span> <span class="special">);</span>
</pre>
<p>
      </p>
<p>
        Note, that we have used <a href="http://www.boost.org/doc/libs/release/libs/phoenix/" target="_top">Boost.Phoenix</a>,
        a great functional programming library, to create and compose the observer.
      </p>
<p>
        The full example can be found here: <a href="https://github.com/headmyshoulder/odeint-v2/blob/master/examples/stiff_system.cpp" target="_top">stiff_system.cpp</a>
      </p>
</div>
<table xmlns:rev="http://www.cs.rpi.edu/~gregod/boost/tools/doc/revision" width="100%"><tr>
<td align="left"></td>
<td align="right"><div class="copyright-footer">Copyright &#169; 2009-2012 Karsten
      Ahnert and Mario Mulansky<p>
        Distributed under the Boost Software License, Version 1.0. (See accompanying
        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>)
      </p>
</div></td>
</tr></table>
<hr>
<div class="spirit-nav">
<a accesskey="p" href="chaotic_systems_and_lyapunov_exponents.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="complex_state_types.html"><img src="../../../../../../../doc/src/images/next.png" alt="Next"></a>
</div>
</body>
</html>