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 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403
|
<!DOCTYPE html>
<html>
<!-- Created by GNU Texinfo 7.1.1, https://www.gnu.org/software/texinfo/ -->
<head>
<meta http-equiv="Content-Type" content="text/html; charset=utf-8">
<title>Solvers (GNU Octave (version 10.3.0))</title>
<meta name="description" content="Solvers (GNU Octave (version 10.3.0))">
<meta name="keywords" content="Solvers (GNU Octave (version 10.3.0))">
<meta name="resource-type" content="document">
<meta name="distribution" content="global">
<meta name="Generator" content="makeinfo">
<meta name="viewport" content="width=device-width,initial-scale=1">
<link href="index.html" rel="start" title="Top">
<link href="Concept-Index.html" rel="index" title="Concept Index">
<link href="index.html#SEC_Contents" rel="contents" title="Table of Contents">
<link href="Nonlinear-Equations.html" rel="up" title="Nonlinear Equations">
<link href="Minimizers.html" rel="next" title="Minimizers">
<style type="text/css">
<!--
a.copiable-link {visibility: hidden; text-decoration: none; line-height: 0em}
div.example {margin-left: 3.2em}
span:hover a.copiable-link {visibility: visible}
strong.def-name {font-family: monospace; font-weight: bold; font-size: larger}
ul.mark-bullet {list-style-type: disc}
-->
</style>
<link rel="stylesheet" type="text/css" href="octave.css">
</head>
<body lang="en">
<div class="section-level-extent" id="Solvers">
<div class="nav-panel">
<p>
Next: <a href="Minimizers.html" accesskey="n" rel="next">Minimizers</a>, Up: <a href="Nonlinear-Equations.html" accesskey="u" rel="up">Nonlinear Equations</a> [<a href="index.html#SEC_Contents" title="Table of contents" rel="contents">Contents</a>][<a href="Concept-Index.html" title="Index" rel="index">Index</a>]</p>
</div>
<hr>
<h3 class="section" id="Solvers-1"><span>20.1 Solvers<a class="copiable-link" href="#Solvers-1"> ¶</a></span></h3>
<p>Octave can solve sets of nonlinear equations of the form
</p>
<div class="example">
<pre class="example-preformatted">F (x) = 0
</pre></div>
<p>using the function <code class="code">fsolve</code>, which is based on the <small class="sc">MINPACK</small>
subroutine <code class="code">hybrd</code>. This is an iterative technique so a starting
point must be provided. This also has the consequence that
convergence is not guaranteed even if a solution exists.
</p>
<a class="anchor" id="XREFfsolve"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-fsolve"><span><code class="def-type"><var class="var">x</var> =</code> <strong class="def-name">fsolve</strong> <code class="def-code-arguments">(<var class="var">fcn</var>, <var class="var">x0</var>)</code><a class="copiable-link" href="#index-fsolve"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-fsolve-1"><span><code class="def-type"><var class="var">x</var> =</code> <strong class="def-name">fsolve</strong> <code class="def-code-arguments">(<var class="var">fcn</var>, <var class="var">x0</var>, <var class="var">options</var>)</code><a class="copiable-link" href="#index-fsolve-1"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-fsolve-2"><span><code class="def-type">[<var class="var">x</var>, <var class="var">fval</var>] =</code> <strong class="def-name">fsolve</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-fsolve-2"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-fsolve-3"><span><code class="def-type">[<var class="var">x</var>, <var class="var">fval</var>, <var class="var">info</var>] =</code> <strong class="def-name">fsolve</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-fsolve-3"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-fsolve-4"><span><code class="def-type">[<var class="var">x</var>, <var class="var">fval</var>, <var class="var">info</var>, <var class="var">output</var>] =</code> <strong class="def-name">fsolve</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-fsolve-4"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-fsolve-5"><span><code class="def-type">[<var class="var">x</var>, <var class="var">fval</var>, <var class="var">info</var>, <var class="var">output</var>, <var class="var">fjac</var>] =</code> <strong class="def-name">fsolve</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-fsolve-5"> ¶</a></span></dt>
<dd><p>Solve a system of nonlinear equations defined by the function <var class="var">fcn</var>.
</p>
<p><var class="var">fcn</var> is a function handle, inline function, or string containing the
name of the function to evaluate. <var class="var">fcn</var> should accept a vector (array)
defining the unknown variables, and return a vector of left-hand sides of
the equations. Right-hand sides are defined to be zeros. In other words,
this function attempts to determine a vector <var class="var">x</var> such that
<code class="code"><var class="var">fcn</var> (<var class="var">x</var>)</code> gives (approximately) all zeros.
</p>
<p><var class="var">x0</var> is an initial guess for the solution. The shape of <var class="var">x0</var> is
preserved in all calls to <var class="var">fcn</var>, but otherwise is treated as a column
vector.
</p>
<p><var class="var">options</var> is a structure specifying additional parameters which
control the algorithm. Currently, <code class="code">fsolve</code> recognizes these options:
<code class="code">"AutoScaling"</code>, <code class="code">"ComplexEqn"</code>, <code class="code">"FinDiffType"</code>,
<code class="code">"FunValCheck"</code>, <code class="code">"Jacobian"</code>, <code class="code">"MaxFunEvals"</code>,
<code class="code">"MaxIter"</code>, <code class="code">"OutputFcn"</code>, <code class="code">"TolFun"</code>, <code class="code">"TolX"</code>,
<code class="code">"TypicalX"</code>, and <code class="code">"Updating"</code>.
</p>
<p>If <code class="code">"AutoScaling"</code> is <code class="code">"on"</code>, the variables will be
automatically scaled according to the column norms of the (estimated)
Jacobian. As a result, <code class="code">"TolFun"</code> becomes scaling-independent. By
default, this option is <code class="code">"off"</code> because it may sometimes deliver
unexpected (though mathematically correct) results.
</p>
<p>If <code class="code">"ComplexEqn"</code> is <code class="code">"on"</code>, <code class="code">fsolve</code> will attempt to solve
complex equations in complex variables, assuming that the equations possess
a complex derivative (i.e., are holomorphic). If this is not what you want,
you should unpack the real and imaginary parts of the system to get a real
system.
</p>
<p>If <code class="code">"Jacobian"</code> is <code class="code">"on"</code>, it specifies that <var class="var">fcn</var>—when
called with 2 output arguments—also returns the Jacobian matrix of
right-hand sides at the requested point.
</p>
<p><code class="code">"MaxFunEvals"</code> proscribes the maximum number of function evaluations
before optimization is halted. The default value is
<code class="code">100 * number_of_variables</code>, i.e., <code class="code">100 * length (<var class="var">x0</var>)</code>.
The value must be a positive integer.
</p>
<p>If <code class="code">"Updating"</code> is <code class="code">"on"</code>, the function will attempt to use
Broyden updates to update the Jacobian, in order to reduce the
number of Jacobian calculations. If your user function always calculates
the Jacobian (regardless of number of output arguments) then this option
provides no advantage and should be disabled.
</p>
<p><code class="code">"TolX"</code> specifies the termination tolerance in the unknown variables,
while <code class="code">"TolFun"</code> is a tolerance for equations. Default is <code class="code">1e-6</code>
for both <code class="code">"TolX"</code> and <code class="code">"TolFun"</code>.
</p>
<p>For a description of the other options,
see <a class="pxref" href="Linear-Least-Squares.html#XREFoptimset"><code class="code">optimset</code></a>. To initialize an options structure
with default values for <code class="code">fsolve</code> use
<code class="code">options = optimset ("fsolve")</code>.
</p>
<p>The first output <var class="var">x</var> is the solution while the second output <var class="var">fval</var>
contains the value of the function <var class="var">fcn</var> evaluated at <var class="var">x</var> (ideally
a vector of all zeros).
</p>
<p>The third output <var class="var">info</var> reports whether the algorithm succeeded and may
take one of the following values:
</p>
<dl class="table">
<dt>1</dt>
<dd><p>Converged to a solution point. Relative residual error is less than
specified by <code class="code">TolFun</code>.
</p>
</dd>
<dt>2</dt>
<dd><p>Last relative step size was less than <code class="code">TolX</code>.
</p>
</dd>
<dt>3</dt>
<dd><p>Last relative decrease in residual was less than <code class="code">TolFun</code>.
</p>
</dd>
<dt>0</dt>
<dd><p>Iteration limit (either <code class="code">MaxIter</code> or <code class="code">MaxFunEvals</code>) exceeded.
</p>
</dd>
<dt>-1</dt>
<dd><p>Stopped by <code class="code">OutputFcn</code>.
</p>
</dd>
<dt>-2</dt>
<dd><p>The Jacobian became excessively small and the search stalled.
</p>
</dd>
<dt>-3</dt>
<dd><p>The trust region radius became excessively small.
</p></dd>
</dl>
<p><var class="var">output</var> is a structure containing runtime information about the
<code class="code">fsolve</code> algorithm. Fields in the structure are:
</p>
<dl class="table">
<dt><code class="code">iterations</code></dt>
<dd><p>Number of iterations through loop.
</p>
</dd>
<dt><code class="code">successful</code></dt>
<dd><p>Number of successful iterations.
</p>
</dd>
<dt><code class="code">funcCount</code></dt>
<dd><p>Number of function evaluations.
</p>
</dd>
</dl>
<p>The final output <var class="var">fjac</var> contains the value of the Jacobian evaluated
at <var class="var">x</var>.
</p>
<p>Note: If you only have a single nonlinear equation of one variable, using
<code class="code">fzero</code> is usually a much better idea.
</p>
<p>Note about user-supplied Jacobians:
As an inherent property of the algorithm, a Jacobian is always requested for
a solution vector whose residual vector is already known, and it is the last
accepted successful step. Often this will be one of the last two calls, but
not always. If the savings by reusing intermediate results from residual
calculation in Jacobian calculation are significant, the best strategy is to
employ <code class="code">OutputFcn</code>: After a vector is evaluated for residuals, if
<code class="code">OutputFcn</code> is called with that vector, then the intermediate results
should be saved for future Jacobian evaluation, and should be kept until a
Jacobian evaluation is requested or until <code class="code">OutputFcn</code> is called with a
different vector, in which case they should be dropped in favor of this most
recent vector. A short example how this can be achieved follows:
</p>
<div class="example">
<pre class="example-preformatted">function [fval, fjac] = user_fcn (x, optimvalues, state)
persistent sav = [], sav0 = [];
if (nargin == 1)
## evaluation call
if (nargout == 1)
sav0.x = x; # mark saved vector
## calculate fval, save results to sav0.
elseif (nargout == 2)
## calculate fjac using sav.
endif
else
## outputfcn call.
if (all (x == sav0.x))
sav = sav0;
endif
## maybe output iteration status, etc.
endif
endfunction
## ...
fsolve (@user_fcn, x0, optimset ("OutputFcn", @user_fcn, ...))
</pre></div>
<p><strong class="strong">See also:</strong> <a class="ref" href="#XREFfzero">fzero</a>, <a class="ref" href="Linear-Least-Squares.html#XREFoptimset">optimset</a>.
</p></dd></dl>
<p>The following is a complete example. To solve the set of equations
</p>
<div class="example">
<div class="group"><pre class="example-preformatted">-2x^2 + 3xy + 4 sin(y) = 6
3x^2 - 2xy^2 + 3 cos(x) = -4
</pre></div></div>
<p>you first need to write a function to compute the value of the given
function. For example:
</p>
<div class="example">
<div class="group"><pre class="example-preformatted">function y = f (x)
y = zeros (2, 1);
y(1) = -2*x(1)^2 + 3*x(1)*x(2) + 4*sin(x(2)) - 6;
y(2) = 3*x(1)^2 - 2*x(1)*x(2)^2 + 3*cos(x(1)) + 4;
endfunction
</pre></div></div>
<p>Then, call <code class="code">fsolve</code> with a specified initial condition to find the
roots of the system of equations. For example, given the function
<code class="code">f</code> defined above,
</p>
<div class="example">
<pre class="example-preformatted">[x, fval, info] = fsolve (@f, [1; 2])
</pre></div>
<p>results in the solution
</p>
<div class="example">
<div class="group"><pre class="example-preformatted">x =
0.57983
2.54621
fval =
-5.7184e-10
5.5460e-10
info = 1
</pre></div></div>
<p>A value of <code class="code">info = 1</code> indicates that the solution has converged.
</p>
<p>When no Jacobian is supplied (as in the example above) it is approximated
numerically. This requires more function evaluations, and hence is
less efficient. In the example above we could compute the Jacobian
analytically as
</p>
<div class="example">
<div class="group"><pre class="example-preformatted">function [y, jac] = f (x)
y = zeros (2, 1);
y(1) = -2*x(1)^2 + 3*x(1)*x(2) + 4*sin(x(2)) - 6;
y(2) = 3*x(1)^2 - 2*x(1)*x(2)^2 + 3*cos(x(1)) + 4;
if (nargout == 2)
jac = zeros (2, 2);
jac(1,1) = 3*x(2) - 4*x(1);
jac(1,2) = 4*cos(x(2)) + 3*x(1);
jac(2,1) = -2*x(2)^2 - 3*sin(x(1)) + 6*x(1);
jac(2,2) = -4*x(1)*x(2);
endif
endfunction
</pre></div></div>
<p>The Jacobian can then be used with the following call to <code class="code">fsolve</code>:
</p>
<div class="example">
<pre class="example-preformatted">[x, fval, info] = fsolve (@f, [1; 2], optimset ("jacobian", "on"));
</pre></div>
<p>which gives the same solution as before.
</p>
<a class="anchor" id="XREFfzero"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-fzero"><span><code class="def-type"><var class="var">x</var> =</code> <strong class="def-name">fzero</strong> <code class="def-code-arguments">(<var class="var">fcn</var>, <var class="var">x0</var>)</code><a class="copiable-link" href="#index-fzero"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-fzero-1"><span><code class="def-type"><var class="var">x</var> =</code> <strong class="def-name">fzero</strong> <code class="def-code-arguments">(<var class="var">fcn</var>, <var class="var">x0</var>, <var class="var">options</var>)</code><a class="copiable-link" href="#index-fzero-1"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-fzero-2"><span><code class="def-type">[<var class="var">x</var>, <var class="var">fval</var>] =</code> <strong class="def-name">fzero</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-fzero-2"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-fzero-3"><span><code class="def-type">[<var class="var">x</var>, <var class="var">fval</var>, <var class="var">info</var>] =</code> <strong class="def-name">fzero</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-fzero-3"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-fzero-4"><span><code class="def-type">[<var class="var">x</var>, <var class="var">fval</var>, <var class="var">info</var>, <var class="var">output</var>] =</code> <strong class="def-name">fzero</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-fzero-4"> ¶</a></span></dt>
<dd><p>Find a zero of a univariate function.
</p>
<p><var class="var">fcn</var> is a function handle, inline function, or string containing the
name of the function to evaluate.
</p>
<p><var class="var">x0</var> should be a two-element vector specifying two points which
bracket a zero. In other words, there must be a change in sign of the
function between <var class="var">x0</var>(1) and <var class="var">x0</var>(2). More mathematically, the
following must hold
</p>
<div class="example">
<pre class="example-preformatted">sign (<var class="var">fcn</var>(<var class="var">x0</var>(1))) * sign (<var class="var">fcn</var>(<var class="var">x0</var>(2))) <= 0
</pre></div>
<p>If <var class="var">x0</var> is a single scalar then several nearby and distant values are
probed in an attempt to obtain a valid bracketing. If this is not
successful, the function fails.
</p>
<p><var class="var">options</var> is a structure specifying additional options. Currently,
<code class="code">fzero</code> recognizes these options:
<code class="code">"Display"</code>, <code class="code">"FunValCheck"</code>, <code class="code">"MaxFunEvals"</code>,
<code class="code">"MaxIter"</code>, <code class="code">"OutputFcn"</code>, and <code class="code">"TolX"</code>.
</p>
<p><code class="code">"MaxFunEvals"</code> proscribes the maximum number of function evaluations
before the search is halted. The default value is <code class="code">Inf</code>.
The value must be a positive integer.
</p>
<p><code class="code">"MaxIter"</code> proscribes the maximum number of algorithm iterations
before the search is halted. The default value is <code class="code">Inf</code>.
The value must be a positive integer.
</p>
<p><code class="code">"TolX"</code> specifies the termination tolerance for the solution <var class="var">x</var>.
The default value is <code class="code">eps</code>.
</p>
<p>For a description of the other options,
see <a class="pxref" href="Linear-Least-Squares.html#XREFoptimset"><code class="code">optimset</code></a>.
To initialize an options structure with default values for <code class="code">fzero</code> use
<code class="code">options = optimset ("fzero")</code>.
</p>
<p>On exit, the function returns <var class="var">x</var>, the approximate zero point, and
<var class="var">fval</var>, the function evaluated at <var class="var">x</var>.
</p>
<p>The third output <var class="var">info</var> reports whether the algorithm succeeded and
may take one of the following values:
</p>
<ul class="itemize mark-bullet">
<li>1
The algorithm converged to a solution.
</li><li>0
Maximum number of iterations or function evaluations has been reached.
</li><li>-1
The algorithm has been terminated by a user <code class="code">OutputFcn</code>.
</li><li>-5
The algorithm may have converged to a singular point.
</li></ul>
<p><var class="var">output</var> is a structure containing runtime information about the
<code class="code">fzero</code> algorithm. Fields in the structure are:
</p>
<ul class="itemize mark-bullet">
<li>iterations
Number of iterations through loop.
</li><li>funcCount
Number of function evaluations.
</li><li>algorithm
The string <code class="code">"bisection, interpolation"</code>.
</li><li>bracketx
A two-element vector with the final bracketing of the zero along the
x-axis.
</li><li>brackety
A two-element vector with the final bracketing of the zero along the
y-axis.
</li></ul>
<p><strong class="strong">See also:</strong> <a class="ref" href="Linear-Least-Squares.html#XREFoptimset">optimset</a>, <a class="ref" href="#XREFfsolve">fsolve</a>.
</p></dd></dl>
</div>
<hr>
<div class="nav-panel">
<p>
Next: <a href="Minimizers.html">Minimizers</a>, Up: <a href="Nonlinear-Equations.html">Nonlinear Equations</a> [<a href="index.html#SEC_Contents" title="Table of contents" rel="contents">Contents</a>][<a href="Concept-Index.html" title="Index" rel="index">Index</a>]</p>
</div>
</body>
</html>
|