File: Providing-the-multidimensional-system-of-equations-to-solve.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 (187 lines) | stat: -rw-r--r-- 8,807 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
<html lang="en">
<head>
<title>Providing the multidimensional system of equations to solve - 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="Multidimensional-Root_002dFinding.html" title="Multidimensional Root-Finding">
<link rel="prev" href="Initializing-the-Multidimensional-Solver.html" title="Initializing the Multidimensional Solver">
<link rel="next" href="Iteration-of-the-multidimensional-solver.html" title="Iteration of the multidimensional solver">
<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="Providing-the-multidimensional-system-of-equations-to-solve"></a>
Next:&nbsp;<a rel="next" accesskey="n" href="Iteration-of-the-multidimensional-solver.html">Iteration of the multidimensional solver</a>,
Previous:&nbsp;<a rel="previous" accesskey="p" href="Initializing-the-Multidimensional-Solver.html">Initializing the Multidimensional Solver</a>,
Up:&nbsp;<a rel="up" accesskey="u" href="Multidimensional-Root_002dFinding.html">Multidimensional Root-Finding</a>
<hr>
</div>

<h3 class="section">34.3 Providing the function to solve</h3>

<p><a name="index-multidimensional-root-finding_002c-providing-a-function-to-solve-2289"></a>
You must provide n functions of n variables for the root
finders to operate on.  In order to allow for general parameters the
functions are defined by the following data types:

<div class="defun">
&mdash; Data Type: <b>gsl_multiroot_function</b><var><a name="index-gsl_005fmultiroot_005ffunction-2290"></a></var><br>
<blockquote><p>This data type defines a general system of functions with parameters.

          <dl>
<dt><code>int (* f) (const gsl_vector * </code><var>x</var><code>, void * </code><var>params</var><code>, gsl_vector * </code><var>f</var><code>)</code><dd>this function should store the vector result
<!-- {$f(x,\hbox{\it params})$} -->
f(x,params) in <var>f</var> for argument <var>x</var> and parameters <var>params</var>,
returning an appropriate error code if the function cannot be computed.

          <br><dt><code>size_t n</code><dd>the dimension of the system, i.e. the number of components of the
vectors <var>x</var> and <var>f</var>.

          <br><dt><code>void * params</code><dd>a pointer to the parameters of the function. 
</dl>
        </p></blockquote></div>

<p class="noindent">Here is an example using Powell's test function,

<pre class="example">     f_1(x) = A x_0 x_1 - 1,
     f_2(x) = exp(-x_0) + exp(-x_1) - (1 + 1/A)
</pre>
   <p class="noindent">with A = 10^4.  The following code defines a
<code>gsl_multiroot_function</code> system <code>F</code> which you could pass to a
solver:

<pre class="example">     struct powell_params { double A; };
     
     int
     powell (gsl_vector * x, void * p, gsl_vector * f) {
        struct powell_params * params
          = *(struct powell_params *)p;
        const double A = (params-&gt;A);
        const double x0 = gsl_vector_get(x,0);
        const double x1 = gsl_vector_get(x,1);
     
        gsl_vector_set (f, 0, A * x0 * x1 - 1);
        gsl_vector_set (f, 1, (exp(-x0) + exp(-x1)
                               - (1.0 + 1.0/A)));
        return GSL_SUCCESS
     }
     
     gsl_multiroot_function F;
     struct powell_params params = { 10000.0 };
     
     F.f = &amp;powell;
     F.n = 2;
     F.params = &amp;params;
</pre>
   <div class="defun">
&mdash; Data Type: <b>gsl_multiroot_function_fdf</b><var><a name="index-gsl_005fmultiroot_005ffunction_005ffdf-2291"></a></var><br>
<blockquote><p>This data type defines a general system of functions with parameters and
the corresponding Jacobian matrix of derivatives,

          <dl>
<dt><code>int (* f) (const gsl_vector * </code><var>x</var><code>, void * </code><var>params</var><code>, gsl_vector * </code><var>f</var><code>)</code><dd>this function should store the vector result
<!-- {$f(x,\hbox{\it params})$} -->
f(x,params) in <var>f</var> for argument <var>x</var> and parameters <var>params</var>,
returning an appropriate error code if the function cannot be computed.

          <br><dt><code>int (* df) (const gsl_vector * </code><var>x</var><code>, void * </code><var>params</var><code>, gsl_matrix * </code><var>J</var><code>)</code><dd>this function should store the <var>n</var>-by-<var>n</var> matrix result
<!-- {$J_{ij} = \partial f_i(x,\hbox{\it params}) / \partial x_j$} -->
J_ij = d f_i(x,params) / d x_j in <var>J</var> for argument <var>x</var>
and parameters <var>params</var>, returning an appropriate error code if the
function cannot be computed.

          <br><dt><code>int (* fdf) (const gsl_vector * </code><var>x</var><code>, void * </code><var>params</var><code>, gsl_vector * </code><var>f</var><code>, gsl_matrix * </code><var>J</var><code>)</code><dd>This function should set the values of the <var>f</var> and <var>J</var> as above,
for arguments <var>x</var> and parameters <var>params</var>.  This function
provides an optimization of the separate functions for f(x) and
J(x)&mdash;it is always faster to compute the function and its
derivative at the same time.

          <br><dt><code>size_t n</code><dd>the dimension of the system, i.e. the number of components of the
vectors <var>x</var> and <var>f</var>.

          <br><dt><code>void * params</code><dd>a pointer to the parameters of the function. 
</dl>
        </p></blockquote></div>

<p class="noindent">The example of Powell's test function defined above can be extended to
include analytic derivatives using the following code,

<pre class="example">     int
     powell_df (gsl_vector * x, void * p, gsl_matrix * J)
     {
        struct powell_params * params
          = *(struct powell_params *)p;
        const double A = (params-&gt;A);
        const double x0 = gsl_vector_get(x,0);
        const double x1 = gsl_vector_get(x,1);
        gsl_matrix_set (J, 0, 0, A * x1);
        gsl_matrix_set (J, 0, 1, A * x0);
        gsl_matrix_set (J, 1, 0, -exp(-x0));
        gsl_matrix_set (J, 1, 1, -exp(-x1));
        return GSL_SUCCESS
     }
     
     int
     powell_fdf (gsl_vector * x, void * p,
                 gsl_matrix * f, gsl_matrix * J) {
        struct powell_params * params
          = *(struct powell_params *)p;
        const double A = (params-&gt;A);
        const double x0 = gsl_vector_get(x,0);
        const double x1 = gsl_vector_get(x,1);
     
        const double u0 = exp(-x0);
        const double u1 = exp(-x1);
     
        gsl_vector_set (f, 0, A * x0 * x1 - 1);
        gsl_vector_set (f, 1, u0 + u1 - (1 + 1/A));
     
        gsl_matrix_set (J, 0, 0, A * x1);
        gsl_matrix_set (J, 0, 1, A * x0);
        gsl_matrix_set (J, 1, 0, -u0);
        gsl_matrix_set (J, 1, 1, -u1);
        return GSL_SUCCESS
     }
     
     gsl_multiroot_function_fdf FDF;
     
     FDF.f = &amp;powell_f;
     FDF.df = &amp;powell_df;
     FDF.fdf = &amp;powell_fdf;
     FDF.n = 2;
     FDF.params = 0;
</pre>
   <p class="noindent">Note that the function <code>powell_fdf</code> is able to reuse existing terms
from the function when calculating the Jacobian, thus saving time.

<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>