File: ODE-Example-programs.html

package info (click to toggle)
gsl-ref-html 1.10-1
  • links: PTS
  • area: main
  • in suites: lenny
  • size: 4,496 kB
  • ctags: 2,955
  • sloc: makefile: 33
file content (226 lines) | stat: -rw-r--r-- 8,316 bytes parent folder | download
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
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
<html lang="en">
<head>
<title>ODE Example programs - GNU Scientific Library -- Reference Manual</title>
<meta http-equiv="Content-Type" content="text/html">
<meta name="description" content="GNU Scientific Library -- Reference Manual">
<meta name="generator" content="makeinfo 4.8">
<link title="Top" rel="start" href="index.html#Top">
<link rel="up" href="Ordinary-Differential-Equations.html" title="Ordinary Differential Equations">
<link rel="prev" href="Evolution.html" title="Evolution">
<link rel="next" href="ODE-References-and-Further-Reading.html" title="ODE References and Further Reading">
<link href="http://www.gnu.org/software/texinfo/" rel="generator-home" title="Texinfo Homepage">
<!--
Copyright (C) 1996, 1997, 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007 The GSL Team.

Permission is granted to copy, distribute and/or modify this document
under the terms of the GNU Free Documentation License, Version 1.2 or
any later version published by the Free Software Foundation; with the
Invariant Sections being ``GNU General Public License'' and ``Free Software
Needs Free Documentation'', the Front-Cover text being ``A GNU Manual'',
and with the Back-Cover Text being (a) (see below).  A copy of the
license is included in the section entitled ``GNU Free Documentation
License''.

(a) The Back-Cover Text is: ``You have freedom to copy and modify this
GNU Manual, like GNU software.''-->
<meta http-equiv="Content-Style-Type" content="text/css">
<style type="text/css"><!--
  pre.display { font-family:inherit }
  pre.format  { font-family:inherit }
  pre.smalldisplay { font-family:inherit; font-size:smaller }
  pre.smallformat  { font-family:inherit; font-size:smaller }
  pre.smallexample { font-size:smaller }
  pre.smalllisp    { font-size:smaller }
  span.sc    { font-variant:small-caps }
  span.roman { font-family:serif; font-weight:normal; } 
  span.sansserif { font-family:sans-serif; font-weight:normal; } 
--></style>
</head>
<body>
<div class="node">
<p>
<a name="ODE-Example-programs"></a>
Next:&nbsp;<a rel="next" accesskey="n" href="ODE-References-and-Further-Reading.html">ODE References and Further Reading</a>,
Previous:&nbsp;<a rel="previous" accesskey="p" href="Evolution.html">Evolution</a>,
Up:&nbsp;<a rel="up" accesskey="u" href="Ordinary-Differential-Equations.html">Ordinary Differential Equations</a>
<hr>
</div>

<h3 class="section">25.5 Examples</h3>

<p><a name="index-Van-der-Pol-oscillator_002c-example-2076"></a>The following program solves the second-order nonlinear Van der Pol
oscillator equation,

<pre class="example">     x''(t) + \mu x'(t) (x(t)^2 - 1) + x(t) = 0
</pre>
   <p class="noindent">This can be converted into a first order system suitable for use with
the routines described in this chapter by introducing a separate
variable for the velocity, y = x'(t),

<pre class="example">     x' = y
     y' = -x + \mu y (1-x^2)
</pre>
   <p class="noindent">The program begins by defining functions for these derivatives and
their Jacobian,

<pre class="example"><pre class="verbatim">     #include &lt;stdio.h>
     #include &lt;gsl/gsl_errno.h>
     #include &lt;gsl/gsl_matrix.h>
     #include &lt;gsl/gsl_odeiv.h>
     
     int
     func (double t, const double y[], double f[],
           void *params)
     {
       double mu = *(double *)params;
       f[0] = y[1];
       f[1] = -y[0] - mu*y[1]*(y[0]*y[0] - 1);
       return GSL_SUCCESS;
     }
     
     int
     jac (double t, const double y[], double *dfdy, 
          double dfdt[], void *params)
     {
       double mu = *(double *)params;
       gsl_matrix_view dfdy_mat 
         = gsl_matrix_view_array (dfdy, 2, 2);
       gsl_matrix * m = &amp;dfdy_mat.matrix; 
       gsl_matrix_set (m, 0, 0, 0.0);
       gsl_matrix_set (m, 0, 1, 1.0);
       gsl_matrix_set (m, 1, 0, -2.0*mu*y[0]*y[1] - 1.0);
       gsl_matrix_set (m, 1, 1, -mu*(y[0]*y[0] - 1.0));
       dfdt[0] = 0.0;
       dfdt[1] = 0.0;
       return GSL_SUCCESS;
     }
     
     int
     main (void)
     {
       const gsl_odeiv_step_type * T 
         = gsl_odeiv_step_rk8pd;
     
       gsl_odeiv_step * s 
         = gsl_odeiv_step_alloc (T, 2);
       gsl_odeiv_control * c 
         = gsl_odeiv_control_y_new (1e-6, 0.0);
       gsl_odeiv_evolve * e 
         = gsl_odeiv_evolve_alloc (2);
     
       double mu = 10;
       gsl_odeiv_system sys = {func, jac, 2, &amp;mu};
     
       double t = 0.0, t1 = 100.0;
       double h = 1e-6;
       double y[2] = { 1.0, 0.0 };
     
       while (t &lt; t1)
         {
           int status = gsl_odeiv_evolve_apply (e, c, s,
                                                &amp;sys, 
                                                &amp;t, t1,
                                                &amp;h, y);
     
           if (status != GSL_SUCCESS)
               break;
     
           printf ("%.5e %.5e %.5e\n", t, y[0], y[1]);
         }
     
       gsl_odeiv_evolve_free (e);
       gsl_odeiv_control_free (c);
       gsl_odeiv_step_free (s);
       return 0;
     }
</pre></pre>
   <p class="noindent">For functions with multiple parameters, the appropriate information
can be passed in through the <var>params</var> argument using a pointer to
a struct.

   <p>The main loop of the program evolves the solution from (y, y') =
(1, 0) at t=0 to t=100.  The step-size h is
automatically adjusted by the controller to maintain an absolute
accuracy of <!-- {$10^{-6}$} -->
10^{-6} in the function values <var>y</var>.

<p class="noindent">To obtain the values at user-specified positions, rather than those
chosen automatically by the control function, the main loop can be modified
to advance the solution from one chosen point to the next.  For example, the
following main loop prints the solution at the points t_i = 0,
1, 2, \dots, 100,

<pre class="example">       for (i = 1; i &lt;= 100; i++)
         {
           double ti = i * t1 / 100.0;
     
           while (t &lt; ti)
             {
               gsl_odeiv_evolve_apply (e, c, s,
                                       &amp;sys,
                                       &amp;t, ti, &amp;h,
                                       y);
             }
     
           printf ("%.5e %.5e %.5e\n", t, y[0], y[1]);
         }
</pre>
   <p class="noindent">Note that arbitrary values of t_i can be used for each stage of
the integration.  The equally spaced points in this example are just
used as an illustration.

   <p>It is also possible to work with a non-adaptive integrator, using only
the stepping function itself.  The following program uses the <code>rk4</code>
fourth-order Runge-Kutta stepping function with a fixed stepsize of
0.01,

<pre class="example"><pre class="verbatim">     int
     main (void)
     {
       const gsl_odeiv_step_type * T 
         = gsl_odeiv_step_rk4;
     
       gsl_odeiv_step * s 
         = gsl_odeiv_step_alloc (T, 2);
     
       double mu = 10;
       gsl_odeiv_system sys = {func, jac, 2, &amp;mu};
     
       double t = 0.0, t1 = 100.0;
       double h = 1e-2;
       double y[2] = { 1.0, 0.0 }, y_err[2];
       double dydt_in[2], dydt_out[2];
     
       /* initialise dydt_in from system parameters */
       GSL_ODEIV_FN_EVAL(&amp;sys, t, y, dydt_in);
     
       while (t &lt; t1)
         {
           int status = gsl_odeiv_step_apply (s, t, h, 
                                              y, y_err, 
                                              dydt_in, 
                                              dydt_out, 
                                              &amp;sys);
     
           if (status != GSL_SUCCESS)
               break;
     
           dydt_in[0] = dydt_out[0];
           dydt_in[1] = dydt_out[1];
     
           t += h;
     
           printf ("%.5e %.5e %.5e\n", t, y[0], y[1]);
         }
     
       gsl_odeiv_step_free (s);
       return 0;
     }
</pre></pre>
   <p class="noindent">The derivatives must be initialized for the starting point t=0
before the first step is taken.  Subsequent steps use the output
derivatives <var>dydt_out</var> as inputs to the next step by copying their
values into <var>dydt_in</var>.

<hr>The GNU Scientific Library - a free numerical library licensed under the GNU GPL<br>Back to the <a href="/software/gsl/">GNU Scientific Library Homepage</a></body></html>