annotate scripts/optimization/fsolve.m @ 8562:a6edd5c23cb5

use replacement methods if qrupdate is not available
author Jaroslav Hajek <highegg@gmail.com>
date Thu, 22 Jan 2009 11:10:47 +0100
parents 1cb63ac13934
children c136d313206a
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
1 ## Copyright (C) 2008 VZLU Prague, a.s.
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
2 ##
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
3 ## This file is part of Octave.
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
4 ##
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
5 ## Octave is free software; you can redistribute it and/or modify it
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
6 ## under the terms of the GNU General Public License as published by
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
7 ## the Free Software Foundation; either version 3 of the License, or (at
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
8 ## your option) any later version.
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
9 ##
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
10 ## Octave is distributed in the hope that it will be useful, but
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
11 ## WITHOUT ANY WARRANTY; without even the implied warranty of
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
12 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
13 ## General Public License for more details.
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
14 ##
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
15 ## You should have received a copy of the GNU General Public License
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
16 ## along with Octave; see the file COPYING. If not, see
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
17 ## <http://www.gnu.org/licenses/>.
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
18 ##
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
19 ## Author: Jaroslav Hajek <highegg@gmail.com>
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
20
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
21 ## -*- texinfo -*-
8513
352d3245d4c1 fsolve.m: doc fix
John W. Eaton <jwe@octave.org>
parents: 8507
diff changeset
22 ## @deftypefn{Function File} {} fsolve (@var{fcn}, @var{x0}, @var{options})
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
23 ## @deftypefnx{Function File} {[@var{x}, @var{fvec}, @var{info}, @var{output}, @var{fjac}]} = fsolve (@var{fcn}, @dots{})
8513
352d3245d4c1 fsolve.m: doc fix
John W. Eaton <jwe@octave.org>
parents: 8507
diff changeset
24 ## Solve a system of nonlinear equations defined by the function @var{fcn}.
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
25 ## @var{fcn} should accepts a vector (array) defining the unknown variables,
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
26 ## and return a vector of left-hand sides of the equations. Right-hand sides
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
27 ## are defined to be zeros.
8514
39867b4aca52 fsolve.m: additional doc fix
John W. Eaton <jwe@octave.org>
parents: 8513
diff changeset
28 ## In other words, this function attempts to determine a vector @var{x} such
39867b4aca52 fsolve.m: additional doc fix
John W. Eaton <jwe@octave.org>
parents: 8513
diff changeset
29 ## that @code{@var{fcn} (@var{x})} gives (approximately) all zeros.
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
30 ## @var{x0} determines a starting guess. The shape of @var{x0} is preserved
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
31 ## in all calls to @var{fcn}, but otherwise it is treated as a column vector.
8515
ec2715c76039 fzero.m, fsolve.m: additional doc fixes
John W. Eaton <jwe@octave.org>
parents: 8514
diff changeset
32 ## @var{options} is a structure specifying additional options.
ec2715c76039 fzero.m, fsolve.m: additional doc fixes
John W. Eaton <jwe@octave.org>
parents: 8514
diff changeset
33 ## Currently, @code{fsolve} recognizes these options:
ec2715c76039 fzero.m, fsolve.m: additional doc fixes
John W. Eaton <jwe@octave.org>
parents: 8514
diff changeset
34 ## @code{"FunValCheck"}, @code{"OutputFcn"}, @code{"TolX"},
ec2715c76039 fzero.m, fsolve.m: additional doc fixes
John W. Eaton <jwe@octave.org>
parents: 8514
diff changeset
35 ## @code{"TolFun"}, @code{"MaxIter"}, @code{"MaxFunEvals"}, and
ec2715c76039 fzero.m, fsolve.m: additional doc fixes
John W. Eaton <jwe@octave.org>
parents: 8514
diff changeset
36 ## @code{"Jacobian"}.
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
37 ##
8513
352d3245d4c1 fsolve.m: doc fix
John W. Eaton <jwe@octave.org>
parents: 8507
diff changeset
38 ## If @code{"Jacobian"} is @code{"on"}, it specifies that @var{fcn},
352d3245d4c1 fsolve.m: doc fix
John W. Eaton <jwe@octave.org>
parents: 8507
diff changeset
39 ## called with 2 output arguments, also returns the Jacobian matrix
352d3245d4c1 fsolve.m: doc fix
John W. Eaton <jwe@octave.org>
parents: 8507
diff changeset
40 ## of right-hand sides at the requested point. @code{"TolX"} specifies
352d3245d4c1 fsolve.m: doc fix
John W. Eaton <jwe@octave.org>
parents: 8507
diff changeset
41 ## the termination tolerance in the unknown variables, while
352d3245d4c1 fsolve.m: doc fix
John W. Eaton <jwe@octave.org>
parents: 8507
diff changeset
42 ## @code{"TolFun"} is a tolerance for equations. Default is @code{1e1*eps}
352d3245d4c1 fsolve.m: doc fix
John W. Eaton <jwe@octave.org>
parents: 8507
diff changeset
43 ## for @code{"TolX"} and @code{1e2*eps} for @code{"TolFun"}.
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
44 ## For description of the other options, see @code{optimset}.
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
45 ##
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
46 ## On return, @var{fval} contains the value of the function @var{fcn}
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
47 ## evaluated at @var{x}, and @var{info} may be one of the following values:
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
48 ##
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
49 ## @table @asis
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
50 ## @item 1
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
51 ## Converged to a solution point. Relative residual error is less than specified
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
52 ## by TolFun.
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
53 ## @item 2
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
54 ## Last relative step size was less that TolX.
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
55 ## @item 3
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
56 ## Last relative decrease in residual was less than TolF.
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
57 ## @item 0
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
58 ## Iteration limit exceeded.
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
59 ## @item -3
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
60 ## The trust region radius became excessively small.
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
61 ## @end table
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
62 ##
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
63 ## Note: If you only have a single nonlinear equation of one variable, using
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
64 ## @code{fzero} is usually a much better idea.
8513
352d3245d4c1 fsolve.m: doc fix
John W. Eaton <jwe@octave.org>
parents: 8507
diff changeset
65 ## @seealso{fzero, optimset}
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
66 ## @end deftypefn
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
67
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
68 function [x, fvec, info, output, fjac] = fsolve (fcn, x0, options)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
69
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
70 if (nargin < 3)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
71 options = struct ();
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
72 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
73
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
74 xsiz = size (x0);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
75 n = numel (x0);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
76
8467
77b8d4aa2743 fsolve.m, fzero.m: style fixes; use strcmpi to compare options
John W. Eaton <jwe@octave.org>
parents: 8466
diff changeset
77 has_jac = strcmpi (optimget (options, "Jacobian", "off"), "on");
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
78 maxiter = optimget (options, "MaxIter", Inf);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
79 maxfev = optimget (options, "MaxFunEvals", Inf);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
80 outfcn = optimget (options, "OutputFcn");
8395
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
81 pivoting = optimget (options, "Pivoting", false);
8467
77b8d4aa2743 fsolve.m, fzero.m: style fixes; use strcmpi to compare options
John W. Eaton <jwe@octave.org>
parents: 8466
diff changeset
82 funvalchk = strcmpi (optimget (options, "FunValCheck", "off"), "on");
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
83
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
84 if (funvalchk)
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
85 ## Replace fun with a guarded version.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
86 fun = @(x) guarded_eval (fun, x);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
87 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
88
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
89 ## These defaults are rather stringent. I think that normally, user
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
90 ## prefers accuracy to performance.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
91
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
92 macheps = eps (class (x0));
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
93
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
94 tolx = optimget (options, "TolX", 1e1*macheps);
8467
77b8d4aa2743 fsolve.m, fzero.m: style fixes; use strcmpi to compare options
John W. Eaton <jwe@octave.org>
parents: 8466
diff changeset
95 tolf = optimget (options, "TolFun", 1e2*macheps);
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
96
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
97 factor = 100;
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
98 ## FIXME: TypicalX corresponds to user scaling (???)
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
99 autodg = true;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
100
8467
77b8d4aa2743 fsolve.m, fzero.m: style fixes; use strcmpi to compare options
John W. Eaton <jwe@octave.org>
parents: 8466
diff changeset
101 niter = 1;
77b8d4aa2743 fsolve.m, fzero.m: style fixes; use strcmpi to compare options
John W. Eaton <jwe@octave.org>
parents: 8466
diff changeset
102 nfev = 0;
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
103
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
104 x = x0(:);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
105 info = 0;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
106
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
107 ## Outer loop.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
108 while (niter < maxiter && nfev < maxfev && ! info)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
109
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
110 ## Calculate function value and Jacobian (possibly via FD).
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
111 ## Handle arbitrary shapes of x and f and remember them.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
112 if (has_jac)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
113 [fvec, fjac] = fcn (reshape (x, xsiz));
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
114 nfev ++;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
115 else
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
116 [fvec, fjac] = __fdjac__ (fcn, reshape (x, xsiz));
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
117 nfev += 1 + length (x);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
118 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
119 fsiz = size (fvec);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
120 fvec = fvec(:);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
121 fn = norm (fvec);
8383
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
122 m = length (fvec);
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
123 n = length (x);
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
124 if (m < n)
8468
866492035ecf fsolve.m, fzero.m: undo part of previous change
John W. Eaton <jwe@octave.org>
parents: 8467
diff changeset
125 error ("fsolve:under", "cannot solve underdetermined systems");
8383
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
126 elseif (m > n && niter == 1)
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
127 if (isempty (optimget (options, "TolFun")))
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
128 warning ("an overdetermined system cannot usually be solved exactly; consider specifying the TolFun option");
8383
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
129 endif
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
130 endif
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
131
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
132 ## Get QR factorization of the jacobian.
8395
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
133 if (pivoting)
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
134 [q, r, p] = qr (fjac);
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
135 else
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
136 [q, r] = qr (fjac);
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
137 endif
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
138
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
139 ## Get column norms, use them as scaling factor.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
140 jcn = norm (fjac, 'columns').';
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
141 if (niter == 1)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
142 if (autodg)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
143 dg = jcn;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
144 dg(dg == 0) = 1;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
145 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
146 xn = norm (dg .* x);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
147 delta = factor * xn;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
148 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
149
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
150 ## Rescale if necessary.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
151 if (autodg)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
152 dg = max (dg, jcn);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
153 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
154
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
155 nfail = 0;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
156 nsuc = 0;
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
157 ## Inner loop.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
158 while (niter <= maxiter && nfev < maxfev && ! info)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
159
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
160 qtf = q'*fvec;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
161
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
162 ## Get TR model (dogleg) minimizer
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
163 ## in case of an overdetermined system, get lsq solution.
8395
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
164 s = - __dogleg__ (r(1:n,:), qtf(1:n), dg, delta);
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
165 if (pivoting)
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
166 s = p * s;
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
167 endif
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
168 sn = norm (dg .* s);
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
169
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
170 if (niter == 1)
8395
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
171 delta = min (delta, sn);
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
172 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
173
8395
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
174 fvec1 = fcn (reshape (x + s, xsiz)) (:);
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
175 fn1 = norm (fvec1);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
176
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
177 if (fn1 < fn)
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
178 ## Scaled actual reduction.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
179 actred = 1 - (fn1/fn)^2;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
180 else
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
181 actred = -1;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
182 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
183
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
184 ## Scaled predicted reduction, and ratio.
8395
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
185 if (pivoting)
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
186 w = qtf + r*(p'*s);
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
187 else
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
188 w = qtf + r*s;
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
189 endif
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
190 t = norm (w);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
191 if (t < fn)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
192 prered = 1 - (t/fn)^2;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
193 ratio = actred / prered;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
194 else
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
195 prered = 0;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
196 ratio = 0;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
197 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
198
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
199 ## Update delta.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
200 if (ratio < 0.1)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
201 nsuc = 0;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
202 nfail ++;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
203 delta *= 0.5;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
204 if (delta <= sqrt (macheps)*xn)
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
205 ## Trust region became uselessly small.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
206 info = -3;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
207 break;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
208 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
209 else
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
210 nfail = 0;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
211 nsuc ++;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
212 if (abs (1-ratio) <= 0.1)
8395
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
213 delta = 2*sn;
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
214 elseif (ratio >= 0.5 || nsuc > 1)
8395
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
215 delta = max (delta, 2*sn);
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
216 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
217 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
218
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
219 if (ratio >= 1e-4)
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
220 ## Successful iteration.
8395
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
221 x += s;
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
222 xn = norm (dg .* x);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
223 fvec = fvec1;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
224 fn = fn1;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
225 niter ++;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
226 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
227
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
228 ## Tests for termination conditions. A mysterious place, anything
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
229 ## can happen if you change something here...
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
230
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
231 ## The rule of thumb (which I'm not sure M*b is quite following)
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
232 ## is that for a tolerance that depends on scaling, only 0 makes
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
233 ## sense as a default value. But 0 usually means uselessly long
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
234 ## iterations, so we need scaling-independent tolerances wherever
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
235 ## possible.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
236
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
237 ## FIXME -- why tolf*n*xn? If abs (e) ~ abs(x) * eps is a vector
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
238 ## of perturbations of x, then norm (fjac*e) <= eps*n*xn, i.e. by
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
239 ## tolf ~ eps we demand as much accuracy as we can expect.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
240 if (fn <= tolf*n*xn)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
241 info = 1;
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
242 ## The following tests done only after successful step.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
243 elseif (actred > 0)
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
244 ## This one is classic. Note that we use scaled variables again,
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
245 ## but compare to scaled step, so nothing bad.
8395
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
246 if (sn <= tolx*xn)
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
247 info = 2;
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
248 ## Again a classic one. It seems weird to use the same tolf
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
249 ## for two different tests, but that's what M*b manual appears
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
250 ## to say.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
251 elseif (actred < tolf)
8383
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
252 info = 3;
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
253 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
254 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
255
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
256 ## Criterion for recalculating jacobian.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
257 if (nfail == 2)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
258 break;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
259 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
260
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
261 ## Compute the scaled Broyden update.
8395
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
262 u = (fvec1 - q*w) / sn;
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
263 v = dg .* ((dg .* s) / sn);
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
264 if (pivoting)
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
265 v = p'*v;
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
266 endif
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
267
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
268 ## Update the QR factorization.
8562
a6edd5c23cb5 use replacement methods if qrupdate is not available
Jaroslav Hajek <highegg@gmail.com>
parents: 8550
diff changeset
269 [q, r] = qrupdate (q, r, u, v);
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
270
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
271 endwhile
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
272 endwhile
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
273
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
274 ## Restore original shapes.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
275 x = reshape (x, xsiz);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
276 fvec = reshape (fvec, fsiz);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
277
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
278 output.iterations = niter;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
279 output.funcCount = niter + 2;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
280
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
281 endfunction
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
282
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
283 ## An assistant function that evaluates a function handle and checks for
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
284 ## bad results.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
285 function fx = guarded_eval (fun, x)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
286 fx = fun (x);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
287 if (! all (isreal (fx)))
8468
866492035ecf fsolve.m, fzero.m: undo part of previous change
John W. Eaton <jwe@octave.org>
parents: 8467
diff changeset
288 error ("fsolve:notreal", "fsolve: non-real value encountered");
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
289 elseif (any (isnan (fx)))
8468
866492035ecf fsolve.m, fzero.m: undo part of previous change
John W. Eaton <jwe@octave.org>
parents: 8467
diff changeset
290 error ("fsolve:isnan", "fsolve: NaN value encountered");
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
291 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
292 endfunction
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
293
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
294 %!function retval = f (p)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
295 %! x = p(1);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
296 %! y = p(2);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
297 %! z = p(3);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
298 %! retval = zeros (3, 1);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
299 %! retval(1) = sin(x) + y**2 + log(z) - 7;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
300 %! retval(2) = 3*x + 2**y -z**3 + 1;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
301 %! retval(3) = x + y + z - 5;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
302 %!test
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
303 %! x_opt = [ 0.599054;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
304 %! 2.395931;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
305 %! 2.005014 ];
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
306 %! tol = 1.0e-5;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
307 %! [x, fval, info] = fsolve (@f, [ 0.5; 2.0; 2.5 ]);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
308 %! assert (info > 0);
8383
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
309 %! assert (norm (x - x_opt, Inf) < tol);
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
310 %! assert (norm (fval) < tol);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
311
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
312 %!function retval = f (p)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
313 %! x = p(1);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
314 %! y = p(2);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
315 %! z = p(3);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
316 %! w = p(4);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
317 %! retval = zeros (4, 1);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
318 %! retval(1) = 3*x + 4*y + exp (z + w) - 1.007;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
319 %! retval(2) = 6*x - 4*y + exp (3*z + w) - 11;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
320 %! retval(3) = x^4 - 4*y^2 + 6*z - 8*w - 20;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
321 %! retval(4) = x^2 + 2*y^3 + z - w - 4;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
322 %!test
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
323 %! x_opt = [ -0.767297326653401, 0.590671081117440, 1.47190018629642, -1.52719341133957 ];
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
324 %! tol = 1.0e-5;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
325 %! [x, fval, info] = fsolve (@f, [-1, 1, 2, -1]);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
326 %! assert (info > 0);
8383
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
327 %! assert (norm (x - x_opt, Inf) < tol);
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
328 %! assert (norm (fval) < tol);
8383
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
329
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
330 %!function retval = f (p)
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
331 %! x = p(1);
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
332 %! y = p(2);
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
333 %! z = p(3);
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
334 %! retval = zeros (3, 1);
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
335 %! retval(1) = sin(x) + y**2 + log(z) - 7;
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
336 %! retval(2) = 3*x + 2**y -z**3 + 1;
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
337 %! retval(3) = x + y + z - 5;
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
338 %! retval(4) = x*x + y - z*log(z) - 1.36;
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
339 %!test
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
340 %! x_opt = [ 0.599054;
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
341 %! 2.395931;
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
342 %! 2.005014 ];
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
343 %! tol = 1.0e-5;
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
344 %! [x, fval, info] = fsolve (@f, [ 0.5; 2.0; 2.5 ], optimset ("TolFun", 1e-6));
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
345 %! assert (info > 0);
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
346 %! assert (norm (x - x_opt, Inf) < tol);
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
347 %! assert (norm (fval) < tol);
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
348
8395
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
349 %!function retval = f (p)
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
350 %! x = p(1);
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
351 %! y = p(2);
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
352 %! z = p(3);
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
353 %! retval = zeros (3, 1);
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
354 %! retval(1) = sin(x) + y**2 + log(z) - 7;
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
355 %! retval(2) = 3*x + 2**y -z**3 + 1;
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
356 %! retval(3) = x + y + z - 5;
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
357 %!test
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
358 %! x_opt = [ 0.599054;
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
359 %! 2.395931;
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
360 %! 2.005014 ];
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
361 %! tol = 1.0e-5;
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
362 %! opt = optimset ("Pivoting", true);
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
363 %! [x, fval, info] = fsolve (@f, [ 0.5; 2.0; 2.5 ], opt);
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
364 %! assert (info > 0);
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
365 %! assert (norm (x - x_opt, Inf) < tol);
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
366 %! assert (norm (fval) < tol);