Mercurial > hg > octave-nkf
annotate scripts/optimization/fsolve.m @ 8857:20589a8f1a33
fix fsolve bug
author | Jaroslav Hajek <highegg@gmail.com> |
---|---|
date | Tue, 24 Feb 2009 13:53:14 +0100 |
parents | 5ce12bca4c51 |
children | 22c8272af34b |
rev | line source |
---|---|
8590 | 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 | 21 ## -*- texinfo -*- |
8513 | 22 ## @deftypefn{Function File} {} fsolve (@var{fcn}, @var{x0}, @var{options}) |
8466 | 23 ## @deftypefnx{Function File} {[@var{x}, @var{fvec}, @var{info}, @var{output}, @var{fjac}]} = fsolve (@var{fcn}, @dots{}) |
8513 | 24 ## Solve a system of nonlinear equations defined by the function @var{fcn}. |
8466 | 25 ## @var{fcn} should accepts a vector (array) defining the unknown variables, |
26 ## and return a vector of left-hand sides of the equations. Right-hand sides | |
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 | 30 ## @var{x0} determines a starting guess. The shape of @var{x0} is preserved |
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 | 35 ## @code{"TolFun"}, @code{"MaxIter"}, @code{"MaxFunEvals"}, |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
36 ## @code{"Jacobian"}, @code{"Updating"} and @code{"ComplexEqn"}. |
8466 | 37 ## |
8513 | 38 ## If @code{"Jacobian"} is @code{"on"}, it specifies that @var{fcn}, |
39 ## called with 2 output arguments, also returns the Jacobian matrix | |
40 ## of right-hand sides at the requested point. @code{"TolX"} specifies | |
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 | 45 ## updates to update the Jacobian, in order to reduce the amount of jacobian |
46 ## calculations. | |
47 ## If your user function always calculates the Jacobian (regardless of number | |
48 ## of output arguments), this option provides no advantage and should be set to | |
49 ## false. | |
50 ## | |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
51 ## @code{"ComplexEqn"} is @code{"on"}, @code{fsolve} will attempt to solve |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
52 ## complex equations in complex variables, assuming that the equations posess a |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
53 ## complex derivative (i.e. are holomorphic). If this is not what you want, |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
54 ## should unpack the real and imaginary parts of the system to get a real |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
55 ## system. |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
56 ## |
8466 | 57 ## For description of the other options, see @code{optimset}. |
58 ## | |
59 ## On return, @var{fval} contains the value of the function @var{fcn} | |
60 ## evaluated at @var{x}, and @var{info} may be one of the following values: | |
61 ## | |
62 ## @table @asis | |
63 ## @item 1 | |
64 ## Converged to a solution point. Relative residual error is less than specified | |
65 ## by TolFun. | |
66 ## @item 2 | |
67 ## Last relative step size was less that TolX. | |
68 ## @item 3 | |
69 ## Last relative decrease in residual was less than TolF. | |
70 ## @item 0 | |
71 ## Iteration limit exceeded. | |
72 ## @item -3 | |
73 ## The trust region radius became excessively small. | |
74 ## @end table | |
75 ## | |
76 ## Note: If you only have a single nonlinear equation of one variable, using | |
77 ## @code{fzero} is usually a much better idea. | |
8513 | 78 ## @seealso{fzero, optimset} |
8466 | 79 ## @end deftypefn |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
80 |
8648
ff61b53eb294
optimization: use PKG_ADD: comments instead of PKG_ADD file
John W. Eaton <jwe@octave.org>
parents:
8647
diff
changeset
|
81 ## PKG_ADD: __all_opts__ ("fsolve"); |
ff61b53eb294
optimization: use PKG_ADD: comments instead of PKG_ADD file
John W. Eaton <jwe@octave.org>
parents:
8647
diff
changeset
|
82 |
8596
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
83 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
|
84 |
8596
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
85 ## Get default options if requested. |
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
86 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
|
87 x = optimset ("MaxIter", 400, "MaxFunEvals", Inf, \ |
8647
06f5dd901f30
implement registering of optimization options
Jaroslav Hajek <highegg@gmail.com>
parents:
8616
diff
changeset
|
88 "Jacobian", "off", "TolX", 1.5e-8, "TolFun", 1.5e-8, |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
89 "OutputFcn", [], "Updating", "on", "FunValCheck", "off", |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
90 "ComplexEqn", "off"); |
8596
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
91 return; |
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
92 endif |
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
93 |
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
94 if (nargin < 2 || nargin > 3 || ! ismatrix (x0)) |
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
95 print_usage (); |
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
96 endif |
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
97 |
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
98 if (ischar (fcn)) |
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
99 fcn = str2func (fcn); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
100 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
101 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
102 xsiz = size (x0); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
103 n = numel (x0); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
104 |
8467
77b8d4aa2743
fsolve.m, fzero.m: style fixes; use strcmpi to compare options
John W. Eaton <jwe@octave.org>
parents:
8466
diff
changeset
|
105 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
|
106 maxiter = optimget (options, "MaxIter", 400); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
107 maxfev = optimget (options, "MaxFunEvals", Inf); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
108 outfcn = optimget (options, "OutputFcn"); |
8596
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
109 updating = strcmpi (optimget (options, "Updating", "on"), "on"); |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
110 complexeqn = strcmpi (optimget (options, "ComplexEqn", "off"), "on"); |
8590 | 111 |
8467
77b8d4aa2743
fsolve.m, fzero.m: style fixes; use strcmpi to compare options
John W. Eaton <jwe@octave.org>
parents:
8466
diff
changeset
|
112 funvalchk = strcmpi (optimget (options, "FunValCheck", "off"), "on"); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
113 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
114 if (funvalchk) |
8604 | 115 ## Replace fcn with a guarded version. |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
116 fcn = @(x) guarded_eval (fcn, x, complexeqn); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
117 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
118 |
8466 | 119 ## These defaults are rather stringent. I think that normally, user |
120 ## prefers accuracy to performance. | |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
121 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
122 macheps = eps (class (x0)); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
123 |
8613
38482007c834
relax scaling in fsolve and use class-dependent default tolerances
Jaroslav Hajek <highegg@gmail.com>
parents:
8604
diff
changeset
|
124 tolx = optimget (options, "TolX", sqrt (macheps)); |
38482007c834
relax scaling in fsolve and use class-dependent default tolerances
Jaroslav Hajek <highegg@gmail.com>
parents:
8604
diff
changeset
|
125 tolf = optimget (options, "TolFun", sqrt (macheps)); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
126 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
127 factor = 100; |
8466 | 128 ## FIXME: TypicalX corresponds to user scaling (???) |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
129 autodg = true; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
130 |
8467
77b8d4aa2743
fsolve.m, fzero.m: style fixes; use strcmpi to compare options
John W. Eaton <jwe@octave.org>
parents:
8466
diff
changeset
|
131 niter = 1; |
77b8d4aa2743
fsolve.m, fzero.m: style fixes; use strcmpi to compare options
John W. Eaton <jwe@octave.org>
parents:
8466
diff
changeset
|
132 nfev = 0; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
133 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
134 x = x0(:); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
135 info = 0; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
136 |
8507 | 137 ## Outer loop. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
138 while (niter < maxiter && nfev < maxfev && ! info) |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
139 |
8507 | 140 ## Calculate function value and Jacobian (possibly via FD). |
141 ## Handle arbitrary shapes of x and f and remember them. | |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
142 if (has_jac) |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
143 [fvec, fjac] = fcn (reshape (x, xsiz)); |
8592
dacfd030633a
handle sparse jacobians in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8590
diff
changeset
|
144 ## If the jacobian is sparse, disable Broyden updating. |
dacfd030633a
handle sparse jacobians in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8590
diff
changeset
|
145 if (issparse (fjac)) |
dacfd030633a
handle sparse jacobians in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8590
diff
changeset
|
146 updating = false; |
dacfd030633a
handle sparse jacobians in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8590
diff
changeset
|
147 endif |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
148 nfev ++; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
149 else |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
150 [fvec, fjac] = __fdjac__ (fcn, reshape (x, xsiz)); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
151 nfev += 1 + length (x); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
152 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
153 fsiz = size (fvec); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
154 fvec = fvec(:); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
155 fn = norm (fvec); |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
156 m = length (fvec); |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
157 n = length (x); |
8590 | 158 |
8763
5ce12bca4c51
update comments in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8693
diff
changeset
|
159 ## For square and overdetermined systems, we update a QR |
8590 | 160 ## factorization of the jacobian to avoid solving a full system in each |
161 ## step. In this case, we pass a triangular matrix to __dogleg__. | |
162 useqr = updating && m >= n && n > 10; | |
163 | |
164 if (useqr) | |
8616
3d75d717cbe0
do not pivot by default in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8615
diff
changeset
|
165 ## FIXME: Currently, pivoting is mostly useless because the \ operator |
3d75d717cbe0
do not pivot by default in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8615
diff
changeset
|
166 ## cannot exploit the resulting props of the triangular factor. |
3d75d717cbe0
do not pivot by default in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8615
diff
changeset
|
167 ## Unpivoted QR is significantly faster so it doesn't seem right to pivot |
3d75d717cbe0
do not pivot by default in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8615
diff
changeset
|
168 ## just to get invariance. Original MINPACK didn't pivot either, at least |
3d75d717cbe0
do not pivot by default in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8615
diff
changeset
|
169 ## when qr updating was used. |
3d75d717cbe0
do not pivot by default in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8615
diff
changeset
|
170 [q, r] = qr (fjac, 0); |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
171 endif |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
172 |
8590 | 173 ## Get column norms, use them as scaling factors. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
174 jcn = norm (fjac, 'columns').'; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
175 if (niter == 1) |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
176 if (autodg) |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
177 dg = jcn; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
178 dg(dg == 0) = 1; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
179 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
180 xn = norm (dg .* x); |
8590 | 181 ## FIXME: something better? |
182 delta = max (factor * xn, 1); | |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
183 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
184 |
8507 | 185 ## Rescale if necessary. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
186 if (autodg) |
8613
38482007c834
relax scaling in fsolve and use class-dependent default tolerances
Jaroslav Hajek <highegg@gmail.com>
parents:
8604
diff
changeset
|
187 ## FIXME: the original minpack used the following rescaling strategy: |
38482007c834
relax scaling in fsolve and use class-dependent default tolerances
Jaroslav Hajek <highegg@gmail.com>
parents:
8604
diff
changeset
|
188 ## dg = max (dg, jcn); |
38482007c834
relax scaling in fsolve and use class-dependent default tolerances
Jaroslav Hajek <highegg@gmail.com>
parents:
8604
diff
changeset
|
189 ## but it seems not good if we start with a bad guess yielding jacobian |
38482007c834
relax scaling in fsolve and use class-dependent default tolerances
Jaroslav Hajek <highegg@gmail.com>
parents:
8604
diff
changeset
|
190 ## columns with large norms that later decrease, because the corresponding |
38482007c834
relax scaling in fsolve and use class-dependent default tolerances
Jaroslav Hajek <highegg@gmail.com>
parents:
8604
diff
changeset
|
191 ## variable will still be overscaled. So instead, we only give the old |
38482007c834
relax scaling in fsolve and use class-dependent default tolerances
Jaroslav Hajek <highegg@gmail.com>
parents:
8604
diff
changeset
|
192 ## scaling a small momentum, but do not honor it. |
38482007c834
relax scaling in fsolve and use class-dependent default tolerances
Jaroslav Hajek <highegg@gmail.com>
parents:
8604
diff
changeset
|
193 |
38482007c834
relax scaling in fsolve and use class-dependent default tolerances
Jaroslav Hajek <highegg@gmail.com>
parents:
8604
diff
changeset
|
194 dg = max (0.1*dg, jcn); |
38482007c834
relax scaling in fsolve and use class-dependent default tolerances
Jaroslav Hajek <highegg@gmail.com>
parents:
8604
diff
changeset
|
195 |
38482007c834
relax scaling in fsolve and use class-dependent default tolerances
Jaroslav Hajek <highegg@gmail.com>
parents:
8604
diff
changeset
|
196 ## It also seems that in the case of fast (and inhomogeneously) changing |
38482007c834
relax scaling in fsolve and use class-dependent default tolerances
Jaroslav Hajek <highegg@gmail.com>
parents:
8604
diff
changeset
|
197 ## jacobian, the Broyden updates are of little use, so maybe we could |
38482007c834
relax scaling in fsolve and use class-dependent default tolerances
Jaroslav Hajek <highegg@gmail.com>
parents:
8604
diff
changeset
|
198 ## skip them if a big disproportional change is expected. The question is, |
38482007c834
relax scaling in fsolve and use class-dependent default tolerances
Jaroslav Hajek <highegg@gmail.com>
parents:
8604
diff
changeset
|
199 ## of course, how to define the above terms :) |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
200 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
201 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
202 nfail = 0; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
203 nsuc = 0; |
8507 | 204 ## Inner loop. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
205 while (niter <= maxiter && nfev < maxfev && ! info) |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
206 |
8590 | 207 if (useqr) |
208 tr_mat = r; | |
209 tr_vec = q'*fvec; | |
210 else | |
211 tr_mat = fjac; | |
212 tr_vec = fvec; | |
213 endif | |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
214 |
8590 | 215 ## Get trust-region model (dogleg) minimizer. |
216 if (useqr) | |
217 qtf = q'*fvec; | |
218 s = - __dogleg__ (r, qtf, dg, delta); | |
219 w = qtf + r * s; | |
220 else | |
221 s = - __dogleg__ (fjac, fvec, dg, delta); | |
222 w = fvec + fjac * s; | |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
223 endif |
8590 | 224 |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
225 sn = norm (dg .* s); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
226 if (niter == 1) |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
227 delta = min (delta, sn); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
228 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
229 |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
230 fvec1 = fcn (reshape (x + s, xsiz)) (:); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
231 fn1 = norm (fvec1); |
8590 | 232 nfev ++; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
233 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
234 if (fn1 < fn) |
8507 | 235 ## Scaled actual reduction. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
236 actred = 1 - (fn1/fn)^2; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
237 else |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
238 actred = -1; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
239 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
240 |
8507 | 241 ## Scaled predicted reduction, and ratio. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
242 t = norm (w); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
243 if (t < fn) |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
244 prered = 1 - (t/fn)^2; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
245 ratio = actred / prered; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
246 else |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
247 prered = 0; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
248 ratio = 0; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
249 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
250 |
8507 | 251 ## Update delta. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
252 if (ratio < 0.1) |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
253 nsuc = 0; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
254 nfail ++; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
255 delta *= 0.5; |
8590 | 256 if (delta <= 1e1*macheps*xn) |
8507 | 257 ## Trust region became uselessly small. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
258 info = -3; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
259 break; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
260 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
261 else |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
262 nfail = 0; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
263 nsuc ++; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
264 if (abs (1-ratio) <= 0.1) |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
265 delta = 2*sn; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
266 elseif (ratio >= 0.5 || nsuc > 1) |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
267 delta = max (delta, 2*sn); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
268 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
269 endif |
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 if (ratio >= 1e-4) |
8507 | 272 ## Successful iteration. |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
273 x += s; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
274 xn = norm (dg .* x); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
275 fvec = fvec1; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
276 fn = fn1; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
277 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
278 |
8857 | 279 niter ++; |
280 | |
8590 | 281 ## FIXME: should outputfcn be only called after a successful iteration? |
8615 | 282 if (! isempty (outfcn)) |
8590 | 283 optimvalues.iter = niter; |
284 optimvalues.funccount = nfev; | |
285 optimvalues.fval = fn; | |
286 optimvalues.searchdirection = s; | |
287 state = 'iter'; | |
288 stop = outfcn (x, optimvalues, state); | |
289 if (stop) | |
290 info = -1; | |
291 break; | |
292 endif | |
293 endif | |
294 | |
8466 | 295 ## Tests for termination conditions. A mysterious place, anything |
296 ## can happen if you change something here... | |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
297 |
8466 | 298 ## The rule of thumb (which I'm not sure M*b is quite following) |
299 ## is that for a tolerance that depends on scaling, only 0 makes | |
300 ## sense as a default value. But 0 usually means uselessly long | |
301 ## iterations, so we need scaling-independent tolerances wherever | |
302 ## possible. | |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
303 |
8466 | 304 ## FIXME -- why tolf*n*xn? If abs (e) ~ abs(x) * eps is a vector |
305 ## of perturbations of x, then norm (fjac*e) <= eps*n*xn, i.e. by | |
306 ## 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
|
307 if (fn <= tolf*n*xn) |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
308 info = 1; |
8466 | 309 ## The following tests done only after successful step. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
310 elseif (actred > 0) |
8466 | 311 ## This one is classic. Note that we use scaled variables again, |
312 ## but compare to scaled step, so nothing bad. | |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
313 if (sn <= tolx*xn) |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
314 info = 2; |
8466 | 315 ## Again a classic one. It seems weird to use the same tolf |
316 ## for two different tests, but that's what M*b manual appears | |
317 ## to say. | |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
318 elseif (actred < tolf) |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
319 info = 3; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
320 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
321 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
322 |
8507 | 323 ## Criterion for recalculating jacobian. |
8590 | 324 if (! updating || nfail == 2) |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
325 break; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
326 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
327 |
8507 | 328 ## Compute the scaled Broyden update. |
8590 | 329 if (useqr) |
330 u = (fvec1 - q*w) / sn; | |
331 v = dg .* ((dg .* s) / sn); | |
332 | |
333 ## Update the QR factorization. | |
334 [q, r] = qrupdate (q, r, u, v); | |
335 else | |
336 u = (fvec1 - w); | |
337 v = dg .* ((dg .* s) / sn); | |
338 | |
339 ## update the jacobian | |
340 fjac += u * v'; | |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
341 endif |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
342 endwhile |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
343 endwhile |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
344 |
8507 | 345 ## Restore original shapes. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
346 x = reshape (x, xsiz); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
347 fvec = reshape (fvec, fsiz); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
348 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
349 output.iterations = niter; |
8590 | 350 output.funcCount = nfev; |
8306
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 endfunction |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
353 |
8466 | 354 ## An assistant function that evaluates a function handle and checks for |
355 ## bad results. | |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
356 function [fx, jx] = guarded_eval (fun, x, complexeqn) |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
357 if (nargout > 1) |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
358 [fx, jx] = fun (x); |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
359 else |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
360 fx = fun (x); |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
361 jx = [] |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
362 endif |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
363 |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
364 if (! complexeqn && ! (all (isreal (fx(:))) && all (isreal (jx(:))))) |
8468
866492035ecf
fsolve.m, fzero.m: undo part of previous change
John W. Eaton <jwe@octave.org>
parents:
8467
diff
changeset
|
365 error ("fsolve:notreal", "fsolve: non-real value encountered"); |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
366 elseif (complexeqn && ! (all (isnumeric (fx(:))) && all (isnumeric(jx(:))))) |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
367 error ("fsolve:notnum", "fsolve: non-numeric value encountered"); |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
368 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
|
369 error ("fsolve:isnan", "fsolve: NaN value encountered"); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
370 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
371 endfunction |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
372 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
373 %!function retval = f (p) |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
374 %! x = p(1); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
375 %! y = p(2); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
376 %! z = p(3); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
377 %! retval = zeros (3, 1); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
378 %! retval(1) = sin(x) + y**2 + log(z) - 7; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
379 %! retval(2) = 3*x + 2**y -z**3 + 1; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
380 %! retval(3) = x + y + z - 5; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
381 %!test |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
382 %! x_opt = [ 0.599054; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
383 %! 2.395931; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
384 %! 2.005014 ]; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
385 %! tol = 1.0e-5; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
386 %! [x, fval, info] = fsolve (@f, [ 0.5; 2.0; 2.5 ]); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
387 %! assert (info > 0); |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
388 %! assert (norm (x - x_opt, Inf) < tol); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
389 %! assert (norm (fval) < tol); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
390 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
391 %!function retval = f (p) |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
392 %! x = p(1); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
393 %! y = p(2); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
394 %! z = p(3); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
395 %! w = p(4); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
396 %! retval = zeros (4, 1); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
397 %! retval(1) = 3*x + 4*y + exp (z + w) - 1.007; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
398 %! retval(2) = 6*x - 4*y + exp (3*z + w) - 11; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
399 %! 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
|
400 %! retval(4) = x^2 + 2*y^3 + z - w - 4; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
401 %!test |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
402 %! x_opt = [ -0.767297326653401, 0.590671081117440, 1.47190018629642, -1.52719341133957 ]; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
403 %! tol = 1.0e-5; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
404 %! [x, fval, info] = fsolve (@f, [-1, 1, 2, -1]); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
405 %! assert (info > 0); |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
406 %! assert (norm (x - x_opt, Inf) < tol); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
407 %! assert (norm (fval) < tol); |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
408 |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
409 %!function retval = f (p) |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
410 %! x = p(1); |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
411 %! y = p(2); |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
412 %! z = p(3); |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
413 %! retval = zeros (3, 1); |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
414 %! retval(1) = sin(x) + y**2 + log(z) - 7; |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
415 %! retval(2) = 3*x + 2**y -z**3 + 1; |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
416 %! retval(3) = x + y + z - 5; |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
417 %! 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
|
418 %!test |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
419 %! x_opt = [ 0.599054; |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
420 %! 2.395931; |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
421 %! 2.005014 ]; |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
422 %! tol = 1.0e-5; |
8590 | 423 %! [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
|
424 %! assert (info > 0); |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
425 %! assert (norm (x - x_opt, Inf) < tol); |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
426 %! assert (norm (fval) < tol); |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
427 |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
428 %!function retval = f (p) |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
429 %! x = p(1); |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
430 %! y = p(2); |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
431 %! z = p(3); |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
432 %! retval = zeros (3, 1); |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
433 %! retval(1) = sin(x) + y**2 + log(z) - 7; |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
434 %! retval(2) = 3*x + 2**y -z**3 + 1; |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
435 %! retval(3) = x + y + z - 5; |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
436 %!test |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
437 %! x_opt = [ 0.599054; |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
438 %! 2.395931; |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
439 %! 2.005014 ]; |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
440 %! tol = 1.0e-5; |
8590 | 441 %! opt = optimset ("Updating", "qrp"); |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
442 %! [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
|
443 %! assert (info > 0); |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
444 %! assert (norm (x - x_opt, Inf) < tol); |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
445 %! assert (norm (fval) < tol); |
8590 | 446 |
447 %!test | |
448 %! b0 = 3; | |
449 %! a0 = 0.2; | |
450 %! x = 0:.5:5; | |
451 %! noise = 1e-5 * sin (100*x); | |
452 %! y = exp (-a0*x) + b0 + noise; | |
453 %! c_opt = [a0, b0]; | |
454 %! tol = 1e-5; | |
455 %! | |
456 %! [c, fval, info, output] = fsolve (@(c) (exp(-c(1)*x) + c(2) - y), [0, 0]); | |
457 %! assert (info > 0); | |
458 %! assert (norm (c - c_opt, Inf) < tol); | |
459 %! assert (norm (fval) < norm (noise)); | |
460 | |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
461 |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
462 %!function y = cfun (x) |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
463 %! y(1) = (1+i)*x(1)^2 - (1-i)*x(2) - 2; |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
464 %! y(2) = sqrt (x(1)*x(2)) - (1-2i)*x(3) + (3-4i); |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
465 %! y(3) = x(1) * x(2) - x(3)^2 + (3+2i); |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
466 |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
467 %!test |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
468 %! x_opt = [-1+i, 1-i, 2+i]; |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
469 %! x = [i, 1, 1+i]; |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
470 %! |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
471 %! [x, f, info] = fsolve (@cfun, x, optimset ("ComplexEqn", "on")); |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
472 %! tol = 1e-5; |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
473 %! assert (norm (f) < tol); |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
474 %! assert (norm (x - x_opt, Inf) < tol); |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
475 |