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
|
@c DO NOT EDIT! Generated automatically by munge-texi.
@c Copyright (C) 1996, 1997, 2007, 2008, 2009 John W. Eaton
@c
@c This file is part of Octave.
@c
@c Octave is free software; you can redistribute it and/or modify it
@c under the terms of the GNU General Public License as published by the
@c Free Software Foundation; either version 3 of the License, or (at
@c your option) any later version.
@c
@c Octave is distributed in the hope that it will be useful, but WITHOUT
@c ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
@c FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
@c for more details.
@c
@c You should have received a copy of the GNU General Public License
@c along with Octave; see the file COPYING. If not, see
@c <http://www.gnu.org/licenses/>.
@node Nonlinear Equations
@chapter Nonlinear Equations
@cindex nonlinear equations
@cindex equations, nonlinear
Octave can solve sets of nonlinear equations of the form
@tex
$$
f (x) = 0
$$
@end tex
@ifnottex
@example
F (x) = 0
@end example
@end ifnottex
@noindent
using the function @code{fsolve}, which is based on the @sc{Minpack}
subroutine @code{hybrd}. This is an iterative technique so a starting
point will have to be provided. This also has the consequence that
convergence is not guaranteed even if a solution exists.
@c ./optimization/fsolve.m
@anchor{doc-fsolve}
@deftypefn {Function File} {} fsolve (@var{fcn}, @var{x0}, @var{options})
@deftypefnx {Function File} {[@var{x}, @var{fvec}, @var{info}, @var{output}, @var{fjac}]} = fsolve (@var{fcn}, @dots{})
Solve a system of nonlinear equations defined by the function @var{fcn}.
@var{fcn} should accepts 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{x} such
that @code{@var{fcn} (@var{x})} gives (approximately) all zeros.
@var{x0} determines a starting guess. The shape of @var{x0} is preserved
in all calls to @var{fcn}, but otherwise it is treated as a column vector.
@var{options} is a structure specifying additional options.
Currently, @code{fsolve} recognizes these options:
@code{"FunValCheck"}, @code{"OutputFcn"}, @code{"TolX"},
@code{"TolFun"}, @code{"MaxIter"}, @code{"MaxFunEvals"},
@code{"Jacobian"}, @code{"Updating"} and @code{"ComplexEqn"}.
If @code{"Jacobian"} is @code{"on"}, it specifies that @var{fcn},
called with 2 output arguments, also returns the Jacobian matrix
of right-hand sides at the requested point. @code{"TolX"} specifies
the termination tolerance in the unknown variables, while
@code{"TolFun"} is a tolerance for equations. Default is @code{1e-7}
for both @code{"TolX"} and @code{"TolFun"}.
If @code{"Updating"} is "on", the function will attempt to use Broyden
updates to update the Jacobian, in order to reduce the amount of jacobian
calculations.
If your user function always calculates the Jacobian (regardless of number
of output arguments), this option provides no advantage and should be set to
false.
@code{"ComplexEqn"} is @code{"on"}, @code{fsolve} 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,
should unpack the real and imaginary parts of the system to get a real
system.
For description of the other options, see @code{optimset}.
On return, @var{fval} contains the value of the function @var{fcn}
evaluated at @var{x}, and @var{info} may be one of the following values:
@table @asis
@item 1
Converged to a solution point. Relative residual error is less than specified
by TolFun.
@item 2
Last relative step size was less that TolX.
@item 3
Last relative decrease in residual was less than TolF.
@item 0
Iteration limit exceeded.
@item -3
The trust region radius became excessively small.
@end table
Note: If you only have a single nonlinear equation of one variable, using
@code{fzero} is usually a much better idea.
@seealso{@ref{doc-fzero,,fzero}, @ref{doc-optimset,,optimset}}
Note about user-supplied jacobians:
As an inherent property of the algorithm, 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 OutputFcn: After a vector is evaluated for residuals, if OutputFcn 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 outputfcn 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:
@example
function [fvec, fjac] = user_func (x, optimvalues, state)
persistent sav = [], sav0 = [];
if (nargin == 1)
## evaluation call
if (nargout == 1)
sav0.x = x; # mark saved vector
## calculate fvec, 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
@dots{}.
fsolve (@@user_func, x0, optimset ("OutputFcn", @@user_func, @dots{}))
@end example
@end deftypefn
Here is a complete example. To solve the set of equations
@tex
$$
\eqalign{-2x^2 + 3xy + 4\sin(y) - 6 &= 0\cr
3x^2 - 2xy^2 + 3\cos(x) + 4 &= 0}
$$
@end tex
@ifinfo
@example
@group
-2x^2 + 3xy + 4 sin(y) = 6
3x^2 - 2xy^2 + 3 cos(x) = -4
@end group
@end example
@end ifinfo
@noindent
you first need to write a function to compute the value of the given
function. For example:
@example
@group
function y = f (x)
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
@end group
@end example
Then, call @code{fsolve} with a specified initial condition to find the
roots of the system of equations. For example, given the function
@code{f} defined above,
@example
[x, fval, info] = fsolve (@@f, [1; 2])
@end example
@noindent
results in the solution
@example
@group
x =
0.57983
2.54621
fval =
-5.7184e-10
5.5460e-10
info = 1
@end group
@end example
@noindent
A value of @code{info = 1} indicates that the solution has converged.
The function @code{perror} may be used to print English messages
corresponding to the numeric error codes. For example,
@example
@group
perror ("fsolve", 1)
@print{} solution converged to requested tolerance
@end group
@end example
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
@iftex
@tex
$$
\left[\matrix{ {\partial f_1 \over \partial x_1} &
{\partial f_1 \over \partial x_2} \cr
{\partial f_2 \over \partial x_1} &
{\partial f_2 \over \partial x_2} \cr}\right] =
\left[\matrix{ 3 x_2 - 4 x_1 &
4 \cos(x_2) + 3 x_1 \cr
-2 x_2^2 - 3 \sin(x_1) + 6 x_1 &
-4 x_1 x_2 \cr }\right]
$$
@end tex
and compute it with the following Octave function
@end iftex
@example
@group
function J = jacobian(x)
J(1,1) = 3*x(2) - 4*x(1);
J(1,2) = 4*cos(x(2)) + 3*x(1);
J(2,1) = -2*x(2)^2 - 3*sin(x(1)) + 6*x(1);
J(2,2) = -4*x(1)*x(2);
endfunction
@end group
@end example
@noindent
The Jacobian can then be used with the following call to @code{fsolve}:
@example
[x, fval, info] = fsolve (@{@@f, @@jacobian@}, [1; 2]);
@end example
@noindent
which gives the same solution as before.
@c ./optimization/fzero.m
@anchor{doc-fzero}
@deftypefn {Function File} {[@var{x}, @var{fval}, @var{info}, @var{output}] =} fzero (@var{fun}, @var{x0}, @var{options})
Find a zero point of a univariate function. @var{fun} should be a function
handle or name. @var{x0} specifies a starting point. @var{options} is a
structure specifying additional options. Currently, @code{fzero}
recognizes these options: @code{"FunValCheck"}, @code{"OutputFcn"},
@code{"TolX"}, @code{"MaxIter"}, @code{"MaxFunEvals"}.
For description of these options, see @ref{doc-optimset,,optimset}.
On exit, the function returns @var{x}, the approximate zero point
and @var{fval}, the function value thereof.
@var{info} is an exit flag that can have these values:
@itemize
@item 1
The algorithm converged to a solution.
@item 0
Maximum number of iterations or function evaluations has been exhausted.
@item -1
The algorithm has been terminated from user output function.
@item -2
A general unexpected error.
@item -3
A non-real value encountered.
@item -4
A NaN value encountered.
@end itemize
@seealso{@ref{doc-optimset,,optimset}, @ref{doc-fsolve,,fsolve}}
@end deftypefn
|