File: Stepping-Functions.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 (172 lines) | stat: -rw-r--r-- 9,988 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
<html lang="en">
<head>
<title>Stepping Functions - 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="Defining-the-ODE-System.html" title="Defining the ODE System">
<link rel="next" href="Adaptive-Step_002dsize-Control.html" title="Adaptive Step-size Control">
<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="Stepping-Functions"></a>
Next:&nbsp;<a rel="next" accesskey="n" href="Adaptive-Step_002dsize-Control.html">Adaptive Step-size Control</a>,
Previous:&nbsp;<a rel="previous" accesskey="p" href="Defining-the-ODE-System.html">Defining the ODE System</a>,
Up:&nbsp;<a rel="up" accesskey="u" href="Ordinary-Differential-Equations.html">Ordinary Differential Equations</a>
<hr>
</div>

<h3 class="section">25.2 Stepping Functions</h3>

<p>The lowest level components are the <dfn>stepping functions</dfn> which
advance a solution from time t to t+h for a fixed
step-size h and estimate the resulting local error.

<div class="defun">
&mdash; Function: gsl_odeiv_step * <b>gsl_odeiv_step_alloc</b> (<var>const gsl_odeiv_step_type * T, size_t dim</var>)<var><a name="index-gsl_005fodeiv_005fstep_005falloc-2033"></a></var><br>
<blockquote><p>This function returns a pointer to a newly allocated instance of a
stepping function of type <var>T</var> for a system of <var>dim</var> dimensions. 
</p></blockquote></div>

<div class="defun">
&mdash; Function: int <b>gsl_odeiv_step_reset</b> (<var>gsl_odeiv_step * s</var>)<var><a name="index-gsl_005fodeiv_005fstep_005freset-2034"></a></var><br>
<blockquote><p>This function resets the stepping function <var>s</var>.  It should be used
whenever the next use of <var>s</var> will not be a continuation of a
previous step. 
</p></blockquote></div>

<div class="defun">
&mdash; Function: void <b>gsl_odeiv_step_free</b> (<var>gsl_odeiv_step * s</var>)<var><a name="index-gsl_005fodeiv_005fstep_005ffree-2035"></a></var><br>
<blockquote><p>This function frees all the memory associated with the stepping function
<var>s</var>. 
</p></blockquote></div>

<div class="defun">
&mdash; Function: const char * <b>gsl_odeiv_step_name</b> (<var>const gsl_odeiv_step * s</var>)<var><a name="index-gsl_005fodeiv_005fstep_005fname-2036"></a></var><br>
<blockquote><p>This function returns a pointer to the name of the stepping function. 
For example,

     <pre class="example">          printf ("step method is '%s'\n",
                   gsl_odeiv_step_name (s));
</pre>
        <p class="noindent">would print something like <code>step method is 'rk4'</code>. 
</p></blockquote></div>

<div class="defun">
&mdash; Function: unsigned int <b>gsl_odeiv_step_order</b> (<var>const gsl_odeiv_step * s</var>)<var><a name="index-gsl_005fodeiv_005fstep_005forder-2037"></a></var><br>
<blockquote><p>This function returns the order of the stepping function on the previous
step.  This order can vary if the stepping function itself is adaptive. 
</p></blockquote></div>

<div class="defun">
&mdash; Function: int <b>gsl_odeiv_step_apply</b> (<var>gsl_odeiv_step * s, double t, double h, double y</var>[]<var>, double yerr</var>[]<var>, const double dydt_in</var>[]<var>, double dydt_out</var>[]<var>, const gsl_odeiv_system * dydt</var>)<var><a name="index-gsl_005fodeiv_005fstep_005fapply-2038"></a></var><br>
<blockquote><p>This function applies the stepping function <var>s</var> to the system of
equations defined by <var>dydt</var>, using the step size <var>h</var> to advance
the system from time <var>t</var> and state <var>y</var> to time <var>t</var>+<var>h</var>. 
The new state of the system is stored in <var>y</var> on output, with an
estimate of the absolute error in each component stored in <var>yerr</var>. 
If the argument <var>dydt_in</var> is not null it should point an array
containing the derivatives for the system at time <var>t</var> on input. This
is optional as the derivatives will be computed internally if they are
not provided, but allows the reuse of existing derivative information. 
On output the new derivatives of the system at time <var>t</var>+<var>h</var> will
be stored in <var>dydt_out</var> if it is not null.

        <p>If the user-supplied functions defined in the system <var>dydt</var> return a
status other than <code>GSL_SUCCESS</code> the step will be aborted.  In this
case, the elements of <var>y</var> will be restored to their pre-step values
and the error code from the user-supplied function will be returned. 
The step-size <var>h</var> will be set to the step-size which caused the error. 
If the function is called again with a smaller step-size, e.g. <var>h</var>/10,
it should be possible to get closer to any singularity. 
To distinguish between error codes from the user-supplied functions and
those from <code>gsl_odeiv_step_apply</code> itself, any user-defined return
values should be distinct from the standard GSL error codes. 
</p></blockquote></div>

   <p>The following algorithms are available,

<div class="defun">
&mdash; Step Type: <b>gsl_odeiv_step_rk2</b><var><a name="index-gsl_005fodeiv_005fstep_005frk2-2039"></a></var><br>
<blockquote><p><a name="index-RK2_002c-Runge_002dKutta-method-2040"></a><a name="index-Runge_002dKutta-methods_002c-ordinary-differential-equations-2041"></a>Embedded Runge-Kutta (2, 3) method. 
</p></blockquote></div>

<div class="defun">
&mdash; Step Type: <b>gsl_odeiv_step_rk4</b><var><a name="index-gsl_005fodeiv_005fstep_005frk4-2042"></a></var><br>
<blockquote><p><a name="index-RK4_002c-Runge_002dKutta-method-2043"></a>4th order (classical) Runge-Kutta. 
</p></blockquote></div>

<div class="defun">
&mdash; Step Type: <b>gsl_odeiv_step_rkf45</b><var><a name="index-gsl_005fodeiv_005fstep_005frkf45-2044"></a></var><br>
<blockquote><p><a name="index-Fehlberg-method_002c-differential-equations-2045"></a><a name="index-RKF45_002c-Runge_002dKutta_002dFehlberg-method-2046"></a>Embedded Runge-Kutta-Fehlberg (4, 5) method.  This method is a good
general-purpose integrator. 
</p></blockquote></div>

<div class="defun">
&mdash; Step Type: <b>gsl_odeiv_step_rkck</b><var><a name="index-gsl_005fodeiv_005fstep_005frkck-2047"></a></var><br>
<blockquote><p><a name="index-Runge_002dKutta-Cash_002dKarp-method-2048"></a><a name="index-Cash_002dKarp_002c-Runge_002dKutta-method-2049"></a>Embedded Runge-Kutta Cash-Karp (4, 5) method. 
</p></blockquote></div>

<div class="defun">
&mdash; Step Type: <b>gsl_odeiv_step_rk8pd</b><var><a name="index-gsl_005fodeiv_005fstep_005frk8pd-2050"></a></var><br>
<blockquote><p><a name="index-Runge_002dKutta-Prince_002dDormand-method-2051"></a><a name="index-Prince_002dDormand_002c-Runge_002dKutta-method-2052"></a>Embedded Runge-Kutta Prince-Dormand (8,9) method. 
</p></blockquote></div>

<div class="defun">
&mdash; Step Type: <b>gsl_odeiv_step_rk2imp</b><var><a name="index-gsl_005fodeiv_005fstep_005frk2imp-2053"></a></var><br>
<blockquote><p>Implicit 2nd order Runge-Kutta at Gaussian points. 
</p></blockquote></div>

<div class="defun">
&mdash; Step Type: <b>gsl_odeiv_step_rk4imp</b><var><a name="index-gsl_005fodeiv_005fstep_005frk4imp-2054"></a></var><br>
<blockquote><p>Implicit 4th order Runge-Kutta at Gaussian points. 
</p></blockquote></div>

<div class="defun">
&mdash; Step Type: <b>gsl_odeiv_step_bsimp</b><var><a name="index-gsl_005fodeiv_005fstep_005fbsimp-2055"></a></var><br>
<blockquote><p><a name="index-Bulirsch_002dStoer-method-2056"></a><a name="index-Bader-and-Deuflhard_002c-Bulirsch_002dStoer-method_002e-2057"></a><a name="index-Deuflhard-and-Bader_002c-Bulirsch_002dStoer-method_002e-2058"></a>Implicit Bulirsch-Stoer method of Bader and Deuflhard.  This algorithm
requires the Jacobian. 
</p></blockquote></div>

<div class="defun">
&mdash; Step Type: <b>gsl_odeiv_step_gear1</b><var><a name="index-gsl_005fodeiv_005fstep_005fgear1-2059"></a></var><br>
<blockquote><p><a name="index-Gear-method_002c-differential-equations-2060"></a>M=1 implicit Gear method. 
</p></blockquote></div>

<div class="defun">
&mdash; Step Type: <b>gsl_odeiv_step_gear2</b><var><a name="index-gsl_005fodeiv_005fstep_005fgear2-2061"></a></var><br>
<blockquote><p>M=2 implicit Gear method. 
</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>