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 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247
|
<!DOCTYPE html PUBLIC "-//W3C//DTD HTML 4.01 Transitional//EN" "http://www.w3.org/TR/html4/loose.dtd">
<html>
<!-- Copyright (C) 1996, 1997, 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009, 2010, 2011, 2012, 2013 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.3 or
any later version published by the Free Software Foundation; with no
Invariant Sections and no cover texts. A copy of the license is
included in the section entitled "GNU Free Documentation License". -->
<!-- Created by GNU Texinfo 5.1, http://www.gnu.org/software/texinfo/ -->
<head>
<title>GNU Scientific Library – Reference Manual: ODE Example programs</title>
<meta name="description" content="GNU Scientific Library – Reference Manual: ODE Example programs">
<meta name="keywords" content="GNU Scientific Library – Reference Manual: ODE Example programs">
<meta name="resource-type" content="document">
<meta name="distribution" content="global">
<meta name="Generator" content="makeinfo">
<meta http-equiv="Content-Type" content="text/html; charset=utf-8">
<link href="index.html#Top" rel="start" title="Top">
<link href="Function-Index.html#Function-Index" rel="index" title="Function Index">
<link href="Ordinary-Differential-Equations.html#Ordinary-Differential-Equations" rel="up" title="Ordinary Differential Equations">
<link href="ODE-References-and-Further-Reading.html#ODE-References-and-Further-Reading" rel="next" title="ODE References and Further Reading">
<link href="Driver.html#Driver" rel="previous" title="Driver">
<style type="text/css">
<!--
a.summary-letter {text-decoration: none}
blockquote.smallquotation {font-size: smaller}
div.display {margin-left: 3.2em}
div.example {margin-left: 3.2em}
div.indentedblock {margin-left: 3.2em}
div.lisp {margin-left: 3.2em}
div.smalldisplay {margin-left: 3.2em}
div.smallexample {margin-left: 3.2em}
div.smallindentedblock {margin-left: 3.2em; font-size: smaller}
div.smalllisp {margin-left: 3.2em}
kbd {font-style:oblique}
pre.display {font-family: inherit}
pre.format {font-family: inherit}
pre.menu-comment {font-family: serif}
pre.menu-preformatted {font-family: serif}
pre.smalldisplay {font-family: inherit; font-size: smaller}
pre.smallexample {font-size: smaller}
pre.smallformat {font-family: inherit; font-size: smaller}
pre.smalllisp {font-size: smaller}
span.nocodebreak {white-space:nowrap}
span.nolinebreak {white-space:nowrap}
span.roman {font-family:serif; font-weight:normal}
span.sansserif {font-family:sans-serif; font-weight:normal}
ul.no-bullet {list-style: none}
-->
</style>
</head>
<body lang="en" bgcolor="#FFFFFF" text="#000000" link="#0000FF" vlink="#800080" alink="#FF0000">
<a name="ODE-Example-programs"></a>
<div class="header">
<p>
Next: <a href="ODE-References-and-Further-Reading.html#ODE-References-and-Further-Reading" accesskey="n" rel="next">ODE References and Further Reading</a>, Previous: <a href="Driver.html#Driver" accesskey="p" rel="previous">Driver</a>, Up: <a href="Ordinary-Differential-Equations.html#Ordinary-Differential-Equations" accesskey="u" rel="up">Ordinary Differential Equations</a> [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
<hr>
<a name="Examples-18"></a>
<h3 class="section">26.6 Examples</h3>
<a name="index-Van-der-Pol-oscillator_002c-example"></a>
<p>The following program solves the second-order nonlinear Van der Pol
oscillator equation,
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, <em>v = u'(t)</em>,
The program begins by defining functions for these derivatives and
their Jacobian. The main function uses driver level functions to solve
the problem. The program evolves the solution from <em>(u, v) = (1,
0)</em> at <em>t=0</em> to <em>t=100</em>. The step-size <em>h</em> is
automatically adjusted by the controller to maintain an absolute
accuracy of <em>10^{-6}</em> in the function values <em>(u, v)</em>.
The loop in the example prints the solution at the points
<em>t_i = 1, 2, \dots, 100</em>.
</p>
<div class="example">
<pre class="verbatim">#include <stdio.h>
#include <gsl/gsl_errno.h>
#include <gsl/gsl_matrix.h>
#include <gsl/gsl_odeiv2.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 = &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)
{
double mu = 10;
gsl_odeiv2_system sys = {func, jac, 2, &mu};
gsl_odeiv2_driver * d =
gsl_odeiv2_driver_alloc_y_new (&sys, gsl_odeiv2_step_rk8pd,
1e-6, 1e-6, 0.0);
int i;
double t = 0.0, t1 = 100.0;
double y[2] = { 1.0, 0.0 };
for (i = 1; i <= 100; i++)
{
double ti = i * t1 / 100.0;
int status = gsl_odeiv2_driver_apply (d, &t, ti, y);
if (status != GSL_SUCCESS)
{
printf ("error, return value=%d\n", status);
break;
}
printf ("%.5e %.5e %.5e\n", t, y[0], y[1]);
}
gsl_odeiv2_driver_free (d);
return 0;
}
</pre></div>
<p>The user can work with the lower level functions directly, as in
the following example. In this case an intermediate result is printed
after each successful step instead of equidistant time points.
</p>
<div class="example">
<pre class="verbatim">int
main (void)
{
const gsl_odeiv2_step_type * T
= gsl_odeiv2_step_rk8pd;
gsl_odeiv2_step * s
= gsl_odeiv2_step_alloc (T, 2);
gsl_odeiv2_control * c
= gsl_odeiv2_control_y_new (1e-6, 0.0);
gsl_odeiv2_evolve * e
= gsl_odeiv2_evolve_alloc (2);
double mu = 10;
gsl_odeiv2_system sys = {func, jac, 2, &mu};
double t = 0.0, t1 = 100.0;
double h = 1e-6;
double y[2] = { 1.0, 0.0 };
while (t < t1)
{
int status = gsl_odeiv2_evolve_apply (e, c, s,
&sys,
&t, t1,
&h, y);
if (status != GSL_SUCCESS)
break;
printf ("%.5e %.5e %.5e\n", t, y[0], y[1]);
}
gsl_odeiv2_evolve_free (e);
gsl_odeiv2_control_free (c);
gsl_odeiv2_step_free (s);
return 0;
}
</pre></div>
<p>For functions with multiple parameters, the appropriate information
can be passed in through the <var>params</var> argument in
<code>gsl_odeiv2_system</code> definition (<var>mu</var> in this example) by using
a pointer to a struct.
</p>
<p>It is also possible to work with a non-adaptive integrator, using only
the stepping function itself,
<code>gsl_odeiv2_driver_apply_fixed_step</code> or
<code>gsl_odeiv2_evolve_apply_fixed_step</code>. The following program uses
the driver level function, with fourth-order
Runge-Kutta stepping function with a fixed stepsize of
0.001.
</p>
<div class="example">
<pre class="verbatim">int
main (void)
{
double mu = 10;
gsl_odeiv2_system sys = { func, jac, 2, &mu };
gsl_odeiv2_driver *d =
gsl_odeiv2_driver_alloc_y_new (&sys, gsl_odeiv2_step_rk4,
1e-3, 1e-8, 1e-8);
double t = 0.0;
double y[2] = { 1.0, 0.0 };
int i, s;
for (i = 0; i < 100; i++)
{
s = gsl_odeiv2_driver_apply_fixed_step (d, &t, 1e-3, 1000, y);
if (s != GSL_SUCCESS)
{
printf ("error: driver returned %d\n", s);
break;
}
printf ("%.5e %.5e %.5e\n", t, y[0], y[1]);
}
gsl_odeiv2_driver_free (d);
return s;
}
</pre></div>
<hr>
<div class="header">
<p>
Next: <a href="ODE-References-and-Further-Reading.html#ODE-References-and-Further-Reading" accesskey="n" rel="next">ODE References and Further Reading</a>, Previous: <a href="Driver.html#Driver" accesskey="p" rel="previous">Driver</a>, Up: <a href="Ordinary-Differential-Equations.html#Ordinary-Differential-Equations" accesskey="u" rel="up">Ordinary Differential Equations</a> [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
</body>
</html>
|