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
|
<html lang="en">
<head>
<title>Defining the ODE System - 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="next" href="Stepping-Functions.html" title="Stepping Functions">
<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="Defining-the-ODE-System"></a>
Next: <a rel="next" accesskey="n" href="Stepping-Functions.html">Stepping Functions</a>,
Up: <a rel="up" accesskey="u" href="Ordinary-Differential-Equations.html">Ordinary Differential Equations</a>
<hr>
</div>
<h3 class="section">25.1 Defining the ODE System</h3>
<p>The routines solve the general n-dimensional first-order system,
<pre class="example"> dy_i(t)/dt = f_i(t, y_1(t), ..., y_n(t))
</pre>
<p class="noindent">for i = 1, \dots, n. The stepping functions rely on the vector
of derivatives f_i and the Jacobian matrix,
<!-- {$J_{ij} = \partial f_i(t, y(t)) / \partial y_j$} -->
J_{ij} = df_i(t,y(t)) / dy_j.
A system of equations is defined using the <code>gsl_odeiv_system</code>
datatype.
<div class="defun">
— Data Type: <b>gsl_odeiv_system</b><var><a name="index-gsl_005fodeiv_005fsystem-2032"></a></var><br>
<blockquote><p>This data type defines a general ODE system with arbitrary parameters.
<dl>
<dt><code>int (* function) (double t, const double y[], double dydt[], void * params)</code><dd>This function should store the vector elements
<!-- {$f_i(t,y,\hbox{\it params})$} -->
f_i(t,y,params) in the array <var>dydt</var>,
for arguments (<var>t</var>,<var>y</var>) and parameters <var>params</var>.
The function should return <code>GSL_SUCCESS</code> if the calculation
was completed successfully. Any other return value indicates
an error.
<br><dt><code>int (* jacobian) (double t, const double y[], double * dfdy, double dfdt[], void * params);</code><dd>This function should store the vector of derivative elements
<!-- {$\partial f_i(t,y,params) / \partial t$} -->
df_i(t,y,params)/dt in the array <var>dfdt</var> and the
Jacobian matrix <!-- {$J_{ij}$} -->
J_{ij} in the array
<var>dfdy</var>, regarded as a row-ordered matrix <code>J(i,j) = dfdy[i * dimension + j]</code>
where <code>dimension</code> is the dimension of the system.
The function should return <code>GSL_SUCCESS</code> if the calculation
was completed successfully. Any other return value indicates
an error.
<p>Some of the simpler solver algorithms do not make use of the Jacobian
matrix, so it is not always strictly necessary to provide it (the
<code>jacobian</code> element of the struct can be replaced by a null pointer
for those algorithms). However, it is useful to provide the Jacobian to allow
the solver algorithms to be interchanged—the best algorithms make use
of the Jacobian.
<br><dt><code>size_t dimension;</code><dd>This is the dimension of the system of equations.
<br><dt><code>void * params</code><dd>This is a pointer to the arbitrary parameters of the system.
</dl>
</p></blockquote></div>
<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>
|