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
  
     | 
    
      <!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, 2014, 2015, 2016 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 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 the freedom to copy and modify this
GNU Manual." -->
<!-- Created by GNU Texinfo 5.1, http://www.gnu.org/software/texinfo/ -->
<head>
<title>GNU Scientific Library – Reference Manual: Providing the multidimensional system of equations to solve</title>
<meta name="description" content="GNU Scientific Library – Reference Manual: Providing the multidimensional system of equations to solve">
<meta name="keywords" content="GNU Scientific Library – Reference Manual: Providing the multidimensional system of equations to solve">
<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="Multidimensional-Root_002dFinding.html#Multidimensional-Root_002dFinding" rel="up" title="Multidimensional Root-Finding">
<link href="Iteration-of-the-multidimensional-solver.html#Iteration-of-the-multidimensional-solver" rel="next" title="Iteration of the multidimensional solver">
<link href="Initializing-the-Multidimensional-Solver.html#Initializing-the-Multidimensional-Solver" rel="previous" title="Initializing the Multidimensional Solver">
<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="Providing-the-multidimensional-system-of-equations-to-solve"></a>
<div class="header">
<p>
Next: <a href="Iteration-of-the-multidimensional-solver.html#Iteration-of-the-multidimensional-solver" accesskey="n" rel="next">Iteration of the multidimensional solver</a>, Previous: <a href="Initializing-the-Multidimensional-Solver.html#Initializing-the-Multidimensional-Solver" accesskey="p" rel="previous">Initializing the Multidimensional Solver</a>, Up: <a href="Multidimensional-Root_002dFinding.html#Multidimensional-Root_002dFinding" accesskey="u" rel="up">Multidimensional Root-Finding</a>   [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
<hr>
<a name="Providing-the-function-to-solve-2"></a>
<h3 class="section">36.3 Providing the function to solve</h3>
<a name="index-multidimensional-root-finding_002c-providing-a-function-to-solve"></a>
<p>You must provide <em>n</em> functions of <em>n</em> variables for the root
finders to operate on.  In order to allow for general parameters the
functions are defined by the following data types:
</p>
<dl>
<dt><a name="index-gsl_005fmultiroot_005ffunction"></a>Data Type: <strong>gsl_multiroot_function</strong></dt>
<dd><p>This data type defines a general system of functions with parameters.
</p>
<dl compact="compact">
<dt><code>int (* f) (const gsl_vector * <var>x</var>, void * <var>params</var>, gsl_vector * <var>f</var>)</code></dt>
<dd><p>this function should store the vector result
<em>f(x,params)</em> 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.
</p>
</dd>
<dt><code>size_t n</code></dt>
<dd><p>the dimension of the system, i.e. the number of components of the
vectors <var>x</var> and <var>f</var>.
</p>
</dd>
<dt><code>void * params</code></dt>
<dd><p>a pointer to the parameters of the function.
</p></dd>
</dl>
</dd></dl>
<p>Here is an example using Powell’s test function,
</p>
<div class="example">
<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></div>
<p>with <em>A = 10^4</em>.  The following code defines a
<code>gsl_multiroot_function</code> system <code>F</code> which you could pass to a
solver:
</p>
<div class="example">
<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->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 = &powell;
F.n = 2;
F.params = &params;
</pre></div>
<dl>
<dt><a name="index-gsl_005fmultiroot_005ffunction_005ffdf"></a>Data Type: <strong>gsl_multiroot_function_fdf</strong></dt>
<dd><p>This data type defines a general system of functions with parameters and
the corresponding Jacobian matrix of derivatives,
</p>
<dl compact="compact">
<dt><code>int (* f) (const gsl_vector * <var>x</var>, void * <var>params</var>, gsl_vector * <var>f</var>)</code></dt>
<dd><p>this function should store the vector result
<em>f(x,params)</em> 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.
</p>
</dd>
<dt><code>int (* df) (const gsl_vector * <var>x</var>, void * <var>params</var>, gsl_matrix * <var>J</var>)</code></dt>
<dd><p>this function should store the <var>n</var>-by-<var>n</var> matrix result
<em>J_ij = d f_i(x,params) / d x_j</em> 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.
</p>
</dd>
<dt><code>int (* fdf) (const gsl_vector * <var>x</var>, void * <var>params</var>, gsl_vector * <var>f</var>, gsl_matrix * <var>J</var>)</code></dt>
<dd><p>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 <em>f(x)</em> and
<em>J(x)</em>—it is always faster to compute the function and its
derivative at the same time.
</p>
</dd>
<dt><code>size_t n</code></dt>
<dd><p>the dimension of the system, i.e. the number of components of the
vectors <var>x</var> and <var>f</var>.
</p>
</dd>
<dt><code>void * params</code></dt>
<dd><p>a pointer to the parameters of the function.
</p></dd>
</dl>
</dd></dl>
<p>The example of Powell’s test function defined above can be extended to
include analytic derivatives using the following code,
</p>
<div class="example">
<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->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->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 = &powell_f;
FDF.df = &powell_df;
FDF.fdf = &powell_fdf;
FDF.n = 2;
FDF.params = 0;
</pre></div>
<p>Note that the function <code>powell_fdf</code> is able to reuse existing terms
from the function when calculating the Jacobian, thus saving time.
</p>
<hr>
<div class="header">
<p>
Next: <a href="Iteration-of-the-multidimensional-solver.html#Iteration-of-the-multidimensional-solver" accesskey="n" rel="next">Iteration of the multidimensional solver</a>, Previous: <a href="Initializing-the-Multidimensional-Solver.html#Initializing-the-Multidimensional-Solver" accesskey="p" rel="previous">Initializing the Multidimensional Solver</a>, Up: <a href="Multidimensional-Root_002dFinding.html#Multidimensional-Root_002dFinding" accesskey="u" rel="up">Multidimensional Root-Finding</a>   [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
</body>
</html>
 
     |