<html lang="en"> <head> <title>Solvers - GNU Octave</title> <meta http-equiv="Content-Type" content="text/html"> <meta name="description" content="GNU Octave"> <meta name="generator" content="makeinfo 4.13"> <link title="Top" rel="start" href="index.html#Top"> <link rel="up" href="Nonlinear-Equations.html#Nonlinear-Equations" title="Nonlinear Equations"> <link rel="next" href="Minimizers.html#Minimizers" title="Minimizers"> <link href="http://www.gnu.org/software/texinfo/" rel="generator-home" title="Texinfo Homepage"> <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"> <a name="Solvers"></a> <p> Next: <a rel="next" accesskey="n" href="Minimizers.html#Minimizers">Minimizers</a>, Up: <a rel="up" accesskey="u" href="Nonlinear-Equations.html#Nonlinear-Equations">Nonlinear Equations</a> <hr> </div> <h3 class="section">20.1 Solvers</h3> <p>Octave can solve sets of nonlinear equations of the form <pre class="example"> F (x) = 0 </pre> <p class="noindent">using the function <code>fsolve</code>, which is based on the <span class="sc">minpack</span> subroutine <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. <!-- fsolve scripts/optimization/fsolve.m --> <p><a name="doc_002dfsolve"></a> <div class="defun"> — Function File: <b>fsolve</b> (<var>fcn, x0, options</var>)<var><a name="index-fsolve-2165"></a></var><br> — Function File: [<var>x</var>, <var>fvec</var>, <var>info</var>, <var>output</var>, <var>fjac</var>] = <b>fsolve</b> (<var>fcn, <small class="dots">...</small></var>)<var><a name="index-fsolve-2166"></a></var><br> <blockquote><p>Solve a system of nonlinear equations defined by the function <var>fcn</var>. <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>x</var> such that <var>fcn</var><code> (</code><var>x</var><code>)</code> gives (approximately) all zeros. <var>x0</var> determines a starting guess. The shape of <var>x0</var> is preserved in all calls to <var>fcn</var>, but otherwise it is treated as a column vector. <var>options</var> is a structure specifying additional options. Currently, <code>fsolve</code> recognizes these options: <code>"FunValCheck"</code>, <code>"OutputFcn"</code>, <code>"TolX"</code>, <code>"TolFun"</code>, <code>"MaxIter"</code>, <code>"MaxFunEvals"</code>, <code>"Jacobian"</code>, <code>"Updating"</code>, <code>"ComplexEqn"</code> <code>"TypicalX"</code>, <code>"AutoScaling"</code> and <code>"FinDiffType"</code>. <p>If <code>"Jacobian"</code> is <code>"on"</code>, it specifies that <var>fcn</var>, called with 2 output arguments, also returns the Jacobian matrix of right-hand sides at the requested point. <code>"TolX"</code> specifies the termination tolerance in the unknown variables, while <code>"TolFun"</code> is a tolerance for equations. Default is <code>1e-7</code> for both <code>"TolX"</code> and <code>"TolFun"</code>. <p>If <code>"AutoScaling"</code> is on, the variables will be automatically scaled according to the column norms of the (estimated) Jacobian. As a result, TolF becomes scaling-independent. By default, this option is off, because it may sometimes deliver unexpected (though mathematically correct) results. <p>If <code>"Updating"</code> 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. <p><code>"ComplexEqn"</code> is <code>"on"</code>, <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, should unpack the real and imaginary parts of the system to get a real system. <p>For description of the other options, see <code>optimset</code>. <p>On return, <var>fval</var> contains the value of the function <var>fcn</var> evaluated at <var>x</var>, and <var>info</var> may be one of the following values: <dl> <dt>1<dd>Converged to a solution point. Relative residual error is less than specified by TolFun. <br><dt>2<dd>Last relative step size was less that TolX. <br><dt>3<dd>Last relative decrease in residual was less than TolF. <br><dt>0<dd>Iteration limit exceeded. <br><dt>-3<dd>The trust region radius became excessively small. </dl> <p>Note: If you only have a single nonlinear equation of one variable, using <code>fzero</code> is usually a much better idea. <p>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: <pre class="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 ## ... fsolve (@user_func, x0, optimset ("OutputFcn", @user_func, ...)) </pre> <!-- Texinfo @sp should work but in practice produces ugly results for HTML. --> <!-- A simple blank line produces the correct behavior. --> <!-- @sp 1 --> <p class="noindent"><strong>See also:</strong> <a href="doc_002dfzero.html#doc_002dfzero">fzero</a>, <a href="doc_002doptimset.html#doc_002doptimset">optimset</a>. </p></blockquote></div> <p>The following is a complete example. To solve the set of equations <pre class="example"> -2x^2 + 3xy + 4 sin(y) = 6 3x^2 - 2xy^2 + 3 cos(x) = -4 </pre> <p class="noindent">you first need to write a function to compute the value of the given function. For example: <pre class="example"> 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> <p>Then, call <code>fsolve</code> with a specified initial condition to find the roots of the system of equations. For example, given the function <code>f</code> defined above, <pre class="example"> [x, fval, info] = fsolve (@f, [1; 2]) </pre> <p class="noindent">results in the solution <pre class="example"> x = 0.57983 2.54621 fval = -5.7184e-10 5.5460e-10 info = 1 </pre> <p class="noindent">A value of <code>info = 1</code> indicates that the solution has converged. <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 <pre class="example"> 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> <p class="noindent">The Jacobian can then be used with the following call to <code>fsolve</code>: <pre class="example"> [x, fval, info] = fsolve (@f, [1; 2], optimset ("jacobian", "on")); </pre> <p class="noindent">which gives the same solution as before. <!-- fzero scripts/optimization/fzero.m --> <p><a name="doc_002dfzero"></a> <div class="defun"> — Function File: <b>fzero</b> (<var>fun, x0</var>)<var><a name="index-fzero-2167"></a></var><br> — Function File: <b>fzero</b> (<var>fun, x0, options</var>)<var><a name="index-fzero-2168"></a></var><br> — Function File: [<var>x</var>, <var>fval</var>, <var>info</var>, <var>output</var>] = <b>fzero</b> (<var><small class="dots">...</small></var>)<var><a name="index-fzero-2169"></a></var><br> <blockquote><p>Find a zero of a univariate function. <p><var>fun</var> is a function handle, inline function, or string containing the name of the function to evaluate. <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>x0</var>(1) and <var>x0</var>(2). More mathematically, the following must hold <pre class="example"> sign (<var>fun</var>(<var>x0</var>(1))) * sign (<var>fun</var>(<var>x0</var>(2))) <= 0 </pre> <p>If <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. <var>options</var> is a structure specifying additional options. Currently, <code>fzero</code> recognizes these options: <code>"FunValCheck"</code>, <code>"OutputFcn"</code>, <code>"TolX"</code>, <code>"MaxIter"</code>, <code>"MaxFunEvals"</code>. For a description of these options, see <a href="doc_002doptimset.html#doc_002doptimset">optimset</a>. <p>On exit, the function returns <var>x</var>, the approximate zero point and <var>fval</var>, the function value thereof. <var>info</var> is an exit flag that can have these values: <ul> <li>1 The algorithm converged to a solution. <li>0 Maximum number of iterations or function evaluations has been reached. <li>-1 The algorithm has been terminated from user output function. <li>-5 The algorithm may have converged to a singular point. </ul> <p><var>output</var> is a structure containing runtime information about the <code>fzero</code> algorithm. Fields in the structure are: <ul> <li>iterations Number of iterations through loop. <li>nfev Number of function evaluations. <li>bracketx A two-element vector with the final bracketing of the zero along the x-axis. <li>brackety A two-element vector with the final bracketing of the zero along the y-axis. </ul> <!-- Texinfo @sp should work but in practice produces ugly results for HTML. --> <!-- A simple blank line produces the correct behavior. --> <!-- @sp 1 --> <p class="noindent"><strong>See also:</strong> <a href="doc_002doptimset.html#doc_002doptimset">optimset</a>, <a href="doc_002dfsolve.html#doc_002dfsolve">fsolve</a>. </p></blockquote></div> </body></html>