annotate scripts/optimization/fsolve.m @ 8596:8833c0b18eb2

enable default settings queries in optim funcs
author Jaroslav Hajek <highegg@gmail.com>
date Tue, 27 Jan 2009 08:15:08 +0100
parents dacfd030633a
children 43f831758d42
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
8590
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
1 ## Copyright (C) 2008, 2009 VZLU Prague, a.s.
8306
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"},
8590
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
35 ## @code{"TolFun"}, @code{"MaxIter"}, @code{"MaxFunEvals"},
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
36 ## @code{"Jacobian"} and @code{"Updating"}.
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
8596
8833c0b18eb2 enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents: 8592
diff changeset
42 ## @code{"TolFun"} is a tolerance for equations. Default is @code{1e-7}
8833c0b18eb2 enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents: 8592
diff changeset
43 ## for both @code{"TolX"} and @code{"TolFun"}.
8833c0b18eb2 enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents: 8592
diff changeset
44 ## If @code{"Updating"} is "on", the function will attempt to use Broyden
8590
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
45 ## updates to update the Jacobian, in order to reduce the amount of jacobian
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
46 ## calculations.
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
47 ## If your user function always calculates the Jacobian (regardless of number
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
48 ## of output arguments), this option provides no advantage and should be set to
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
49 ## false.
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
50 ##
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
51 ## For description of the other options, see @code{optimset}.
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
52 ##
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
53 ## 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
54 ## 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
55 ##
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
56 ## @table @asis
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
57 ## @item 1
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
58 ## 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
59 ## by TolFun.
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
60 ## @item 2
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
61 ## Last relative step size was less that TolX.
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
62 ## @item 3
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
63 ## Last relative decrease in residual was less than TolF.
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
64 ## @item 0
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
65 ## Iteration limit exceeded.
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
66 ## @item -3
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
67 ## The trust region radius became excessively small.
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
68 ## @end table
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
69 ##
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
70 ## 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
71 ## @code{fzero} is usually a much better idea.
8513
352d3245d4c1 fsolve.m: doc fix
John W. Eaton <jwe@octave.org>
parents: 8507
diff changeset
72 ## @seealso{fzero, optimset}
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
73 ## @end deftypefn
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
74
8596
8833c0b18eb2 enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents: 8592
diff changeset
75 function [x, fvec, info, output, fjac] = fsolve (fcn, x0, options = struct ())
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
76
8596
8833c0b18eb2 enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents: 8592
diff changeset
77 ## Get default options if requested.
8833c0b18eb2 enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents: 8592
diff changeset
78 if (nargin == 1 && ischar (fcn) && strcmp (fcn, 'defaults'))
8833c0b18eb2 enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents: 8592
diff changeset
79 x = optimset ("MaxIter", 400, "MaxFunEvals", Inf, \
8833c0b18eb2 enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents: 8592
diff changeset
80 "Jacobian", "off", "TolX", 1e-7, "TolF", 1e-7,
8833c0b18eb2 enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents: 8592
diff changeset
81 "OutputFcn", [], "Updating", "on", "FunValCheck", "off");
8833c0b18eb2 enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents: 8592
diff changeset
82 return;
8833c0b18eb2 enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents: 8592
diff changeset
83 endif
8833c0b18eb2 enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents: 8592
diff changeset
84
8833c0b18eb2 enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents: 8592
diff changeset
85 if (nargin < 2 || nargin > 3 || ! ismatrix (x0))
8833c0b18eb2 enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents: 8592
diff changeset
86 print_usage ();
8833c0b18eb2 enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents: 8592
diff changeset
87 endif
8833c0b18eb2 enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents: 8592
diff changeset
88
8833c0b18eb2 enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents: 8592
diff changeset
89 if (ischar (fcn))
8833c0b18eb2 enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents: 8592
diff changeset
90 fcn = str2func (fcn);
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
91 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
92
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
93 xsiz = size (x0);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
94 n = numel (x0);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
95
8467
77b8d4aa2743 fsolve.m, fzero.m: style fixes; use strcmpi to compare options
John W. Eaton <jwe@octave.org>
parents: 8466
diff changeset
96 has_jac = strcmpi (optimget (options, "Jacobian", "off"), "on");
8596
8833c0b18eb2 enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents: 8592
diff changeset
97 maxiter = optimget (options, "MaxIter", 400);
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
98 maxfev = optimget (options, "MaxFunEvals", Inf);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
99 outfcn = optimget (options, "OutputFcn");
8596
8833c0b18eb2 enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents: 8592
diff changeset
100 updating = strcmpi (optimget (options, "Updating", "on"), "on");
8590
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
101
8467
77b8d4aa2743 fsolve.m, fzero.m: style fixes; use strcmpi to compare options
John W. Eaton <jwe@octave.org>
parents: 8466
diff changeset
102 funvalchk = strcmpi (optimget (options, "FunValCheck", "off"), "on");
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 if (funvalchk)
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
105 ## Replace fun with a guarded version.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
106 fun = @(x) guarded_eval (fun, x);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
107 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
108
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
109 ## 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
110 ## prefers accuracy to performance.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
111
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
112 macheps = eps (class (x0));
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
113
8596
8833c0b18eb2 enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents: 8592
diff changeset
114 tolx = optimget (options, "TolX", 1e-7);
8833c0b18eb2 enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents: 8592
diff changeset
115 tolf = optimget (options, "TolFun", 1e-7);
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
116
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
117 factor = 100;
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
118 ## FIXME: TypicalX corresponds to user scaling (???)
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
119 autodg = true;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
120
8467
77b8d4aa2743 fsolve.m, fzero.m: style fixes; use strcmpi to compare options
John W. Eaton <jwe@octave.org>
parents: 8466
diff changeset
121 niter = 1;
77b8d4aa2743 fsolve.m, fzero.m: style fixes; use strcmpi to compare options
John W. Eaton <jwe@octave.org>
parents: 8466
diff changeset
122 nfev = 0;
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
123
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
124 x = x0(:);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
125 info = 0;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
126
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
127 ## Outer loop.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
128 while (niter < maxiter && nfev < maxfev && ! info)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
129
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
130 ## Calculate function value and Jacobian (possibly via FD).
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
131 ## Handle arbitrary shapes of x and f and remember them.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
132 if (has_jac)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
133 [fvec, fjac] = fcn (reshape (x, xsiz));
8592
dacfd030633a handle sparse jacobians in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8590
diff changeset
134 ## If the jacobian is sparse, disable Broyden updating.
dacfd030633a handle sparse jacobians in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8590
diff changeset
135 if (issparse (fjac))
dacfd030633a handle sparse jacobians in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8590
diff changeset
136 updating = false;
dacfd030633a handle sparse jacobians in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8590
diff changeset
137 endif
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
138 nfev ++;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
139 else
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
140 [fvec, fjac] = __fdjac__ (fcn, reshape (x, xsiz));
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
141 nfev += 1 + length (x);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
142 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
143 fsiz = size (fvec);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
144 fvec = fvec(:);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
145 fn = norm (fvec);
8383
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
146 m = length (fvec);
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
147 n = length (x);
8590
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
148
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
149 ## For square and overdetermined systems, we update a (pivoted) QR
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
150 ## factorization of the jacobian to avoid solving a full system in each
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
151 ## step. In this case, we pass a triangular matrix to __dogleg__.
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
152 ## Pivoted QR is used for slightly better robustness and invariance
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
153 ## w.r.t. permutations of variables.
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
154 useqr = updating && m >= n && n > 10;
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
155
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
156 if (useqr)
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
157 ## Get QR factorization of the jacobian, optionally with column pivoting.
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
158 ## TODO: pivoting is only done in the first step, to get invariance
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
159 ## w.r.t. permutations of variables. Maybe it could be beneficial to
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
160 ## repivot occassionally?
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
161 if (niter == 1)
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
162 [q, r, p] = qr (fjac, 0);
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
163 ## p is a column vector. Blame Matlab for the inconsistency.
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
164 ## Octave can handle permutation matrices gracefully.
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
165 p = eye (n)(:, p);
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
166 else
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
167 [q, r] = qr (fjac*p, 0);
8383
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
168 endif
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
169 endif
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
170
8590
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
171 ## Get column norms, use them as scaling factors.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
172 jcn = norm (fjac, 'columns').';
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
173 if (niter == 1)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
174 if (autodg)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
175 dg = jcn;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
176 dg(dg == 0) = 1;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
177 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
178 xn = norm (dg .* x);
8590
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
179 ## FIXME: something better?
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
180 delta = max (factor * xn, 1);
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
181 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
182
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
183 ## Rescale if necessary.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
184 if (autodg)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
185 dg = max (dg, jcn);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
186 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
187
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
188 nfail = 0;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
189 nsuc = 0;
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
190 ## Inner loop.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
191 while (niter <= maxiter && nfev < maxfev && ! info)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
192
8590
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
193 if (useqr)
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
194 tr_mat = r;
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
195 tr_vec = q'*fvec;
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
196 else
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
197 tr_mat = fjac;
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
198 tr_vec = fvec;
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
199 endif
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
200
8590
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
201 ## Get trust-region model (dogleg) minimizer.
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
202 if (useqr)
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
203 qtf = q'*fvec;
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
204 s = - __dogleg__ (r, qtf, dg, delta);
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
205 w = qtf + r * s;
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
206 s = p * s;
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
207 else
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
208 s = - __dogleg__ (fjac, fvec, dg, delta);
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
209 w = fvec + fjac * s;
8395
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
210 endif
8590
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
211
8395
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
212 sn = norm (dg .* s);
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
213 if (niter == 1)
8395
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
214 delta = min (delta, sn);
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
215 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
216
8395
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
217 fvec1 = fcn (reshape (x + s, xsiz)) (:);
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
218 fn1 = norm (fvec1);
8590
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
219 nfev ++;
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
220
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
221 if (fn1 < fn)
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
222 ## Scaled actual reduction.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
223 actred = 1 - (fn1/fn)^2;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
224 else
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
225 actred = -1;
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
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
228 ## Scaled predicted reduction, and ratio.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
229 t = norm (w);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
230 if (t < fn)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
231 prered = 1 - (t/fn)^2;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
232 ratio = actred / prered;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
233 else
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
234 prered = 0;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
235 ratio = 0;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
236 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
237
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
238 ## Update delta.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
239 if (ratio < 0.1)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
240 nsuc = 0;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
241 nfail ++;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
242 delta *= 0.5;
8590
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
243 if (delta <= 1e1*macheps*xn)
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
244 ## Trust region became uselessly small.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
245 info = -3;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
246 break;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
247 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
248 else
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
249 nfail = 0;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
250 nsuc ++;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
251 if (abs (1-ratio) <= 0.1)
8395
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
252 delta = 2*sn;
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
253 elseif (ratio >= 0.5 || nsuc > 1)
8395
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
254 delta = max (delta, 2*sn);
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
255 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
256 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
257
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
258 if (ratio >= 1e-4)
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
259 ## Successful iteration.
8395
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
260 x += s;
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
261 xn = norm (dg .* x);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
262 fvec = fvec1;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
263 fn = fn1;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
264 niter ++;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
265 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
266
8590
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
267 ## FIXME: should outputfcn be only called after a successful iteration?
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
268 if (outfcn)
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
269 optimvalues.iter = niter;
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
270 optimvalues.funccount = nfev;
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
271 optimvalues.fval = fn;
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
272 optimvalues.searchdirection = s;
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
273 state = 'iter';
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
274 stop = outfcn (x, optimvalues, state);
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
275 if (stop)
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
276 info = -1;
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
277 break;
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
278 endif
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
279 endif
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
280
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
281 ## Tests for termination conditions. A mysterious place, anything
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
282 ## can happen if you change something here...
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
283
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
284 ## 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
285 ## 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
286 ## 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
287 ## iterations, so we need scaling-independent tolerances wherever
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
288 ## possible.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
289
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
290 ## 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
291 ## 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
292 ## 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
293 if (fn <= tolf*n*xn)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
294 info = 1;
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
295 ## The following tests done only after successful step.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
296 elseif (actred > 0)
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
297 ## 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
298 ## but compare to scaled step, so nothing bad.
8395
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
299 if (sn <= tolx*xn)
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
300 info = 2;
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
301 ## 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
302 ## 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
303 ## to say.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
304 elseif (actred < tolf)
8383
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
305 info = 3;
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
306 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
307 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
308
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
309 ## Criterion for recalculating jacobian.
8590
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
310 if (! updating || nfail == 2)
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
311 break;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
312 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
313
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
314 ## Compute the scaled Broyden update.
8590
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
315 if (useqr)
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
316 u = (fvec1 - q*w) / sn;
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
317 v = dg .* ((dg .* s) / sn);
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
318 v = p' * v;
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
319
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
320 ## Update the QR factorization.
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
321 [q, r] = qrupdate (q, r, u, v);
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
322 else
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
323 u = (fvec1 - w);
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
324 v = dg .* ((dg .* s) / sn);
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
325
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
326 ## update the jacobian
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
327 fjac += u * v';
8395
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
328 endif
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
329 endwhile
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
330 endwhile
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
331
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8468
diff changeset
332 ## Restore original shapes.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
333 x = reshape (x, xsiz);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
334 fvec = reshape (fvec, fsiz);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
335
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
336 output.iterations = niter;
8590
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
337 output.funcCount = nfev;
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
338
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
339 endfunction
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
340
8466
4d008d9f0ccf fsolve.m: style fixes
John W. Eaton <jwe@octave.org>
parents: 8395
diff changeset
341 ## 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
342 ## bad results.
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
343 function fx = guarded_eval (fun, x)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
344 fx = fun (x);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
345 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
346 error ("fsolve:notreal", "fsolve: non-real value encountered");
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
347 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
348 error ("fsolve:isnan", "fsolve: NaN value encountered");
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
349 endif
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
350 endfunction
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
351
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
352 %!function retval = f (p)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
353 %! x = p(1);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
354 %! y = p(2);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
355 %! z = p(3);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
356 %! retval = zeros (3, 1);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
357 %! retval(1) = sin(x) + y**2 + log(z) - 7;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
358 %! retval(2) = 3*x + 2**y -z**3 + 1;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
359 %! retval(3) = x + y + z - 5;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
360 %!test
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
361 %! x_opt = [ 0.599054;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
362 %! 2.395931;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
363 %! 2.005014 ];
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
364 %! tol = 1.0e-5;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
365 %! [x, fval, info] = fsolve (@f, [ 0.5; 2.0; 2.5 ]);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
366 %! assert (info > 0);
8383
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
367 %! assert (norm (x - x_opt, Inf) < tol);
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
368 %! assert (norm (fval) < tol);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
369
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
370 %!function retval = f (p)
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
371 %! x = p(1);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
372 %! y = p(2);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
373 %! z = p(3);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
374 %! w = p(4);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
375 %! retval = zeros (4, 1);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
376 %! retval(1) = 3*x + 4*y + exp (z + w) - 1.007;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
377 %! retval(2) = 6*x - 4*y + exp (3*z + w) - 11;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
378 %! 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
379 %! retval(4) = x^2 + 2*y^3 + z - w - 4;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
380 %!test
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
381 %! x_opt = [ -0.767297326653401, 0.590671081117440, 1.47190018629642, -1.52719341133957 ];
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
382 %! tol = 1.0e-5;
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
383 %! [x, fval, info] = fsolve (@f, [-1, 1, 2, -1]);
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
384 %! assert (info > 0);
8383
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
385 %! assert (norm (x - x_opt, Inf) < tol);
8306
43795cf108d0 initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
386 %! assert (norm (fval) < tol);
8383
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
387
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
388 %!function retval = f (p)
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
389 %! x = p(1);
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
390 %! y = p(2);
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
391 %! z = p(3);
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
392 %! retval = zeros (3, 1);
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
393 %! retval(1) = sin(x) + y**2 + log(z) - 7;
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
394 %! retval(2) = 3*x + 2**y -z**3 + 1;
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
395 %! retval(3) = x + y + z - 5;
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
396 %! 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
397 %!test
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
398 %! x_opt = [ 0.599054;
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
399 %! 2.395931;
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
400 %! 2.005014 ];
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
401 %! tol = 1.0e-5;
8590
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
402 %! [x, fval, info] = fsolve (@f, [ 0.5; 2.0; 2.5 ]);
8383
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
403 %! assert (info > 0);
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
404 %! assert (norm (x - x_opt, Inf) < tol);
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
405 %! assert (norm (fval) < tol);
a762d9daa700 allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8306
diff changeset
406
8395
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
407 %!function retval = f (p)
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
408 %! x = p(1);
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
409 %! y = p(2);
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
410 %! z = p(3);
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
411 %! retval = zeros (3, 1);
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
412 %! retval(1) = sin(x) + y**2 + log(z) - 7;
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
413 %! retval(2) = 3*x + 2**y -z**3 + 1;
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
414 %! retval(3) = x + y + z - 5;
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
415 %!test
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
416 %! x_opt = [ 0.599054;
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
417 %! 2.395931;
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
418 %! 2.005014 ];
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
419 %! tol = 1.0e-5;
8590
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
420 %! opt = optimset ("Updating", "qrp");
8395
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
421 %! [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
422 %! assert (info > 0);
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
423 %! assert (norm (x - x_opt, Inf) < tol);
88fd356b0d95 optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8383
diff changeset
424 %! assert (norm (fval) < tol);
8590
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
425
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
426 %!test
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
427 %! b0 = 3;
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
428 %! a0 = 0.2;
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
429 %! x = 0:.5:5;
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
430 %! noise = 1e-5 * sin (100*x);
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
431 %! y = exp (-a0*x) + b0 + noise;
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
432 %! c_opt = [a0, b0];
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
433 %! tol = 1e-5;
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
434 %!
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
435 %! [c, fval, info, output] = fsolve (@(c) (exp(-c(1)*x) + c(2) - y), [0, 0]);
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
436 %! assert (info > 0);
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
437 %! assert (norm (c - c_opt, Inf) < tol);
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
438 %! assert (norm (fval) < norm (noise));
c136d313206a improvements to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 8562
diff changeset
439