Mercurial > hg > octave-lyh
annotate scripts/optimization/fsolve.m @ 17392:8c5878260636
doc: Fix typo in fieldnames docstring.
* scripts/general/fieldnames.m: Eliminate stray parenthesis in docstring.
author | Rik <rik@octave.org> |
---|---|
date | Fri, 06 Sep 2013 08:35:59 -0700 |
parents | 1c89599167a6 |
children |
rev | line source |
---|---|
14138
72c96de7a403
maint: update copyright notices for 2012
John W. Eaton <jwe@octave.org>
parents:
13305
diff
changeset
|
1 ## Copyright (C) 2008-2012 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 -*- |
9143
74d5c1a4ca96
Eliminate 'unbalanced parentheses in @def...' error during texi2pdf.
Rik <rdrider0-list@yahoo.com>
parents:
9075
diff
changeset
|
22 ## @deftypefn {Function File} {} fsolve (@var{fcn}, @var{x0}, @var{options}) |
9724
f22bbc5d56e9
Fix various incorrect usages of TeXinfo deffn and deftypefn macros
Rik <rdrider0-list@yahoo.com>
parents:
9628
diff
changeset
|
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}. |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
25 ## @var{fcn} should accept a vector (array) defining the unknown variables, |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8997
diff
changeset
|
26 ## and return a vector of left-hand sides of the equations. Right-hand sides |
8466 | 27 ## are defined to be zeros. |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
28 ## In other words, this function attempts to determine a vector @var{x} such |
8514
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. |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8997
diff
changeset
|
30 ## @var{x0} determines a starting guess. The shape of @var{x0} is preserved |
8466 | 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: |
17289
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17253
diff
changeset
|
34 ## @qcode{"FunValCheck"}, @qcode{"OutputFcn"}, @qcode{"TolX"}, |
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17253
diff
changeset
|
35 ## @qcode{"TolFun"}, @qcode{"MaxIter"}, @qcode{"MaxFunEvals"}, |
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17253
diff
changeset
|
36 ## @qcode{"Jacobian"}, @qcode{"Updating"}, @qcode{"ComplexEqn"} |
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17253
diff
changeset
|
37 ## @qcode{"TypicalX"}, @qcode{"AutoScaling"} and @qcode{"FinDiffType"}. |
8466 | 38 ## |
17289
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17253
diff
changeset
|
39 ## If @qcode{"Jacobian"} is @qcode{"on"}, it specifies that @var{fcn}, |
8513 | 40 ## called with 2 output arguments, also returns the Jacobian matrix |
17289
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17253
diff
changeset
|
41 ## of right-hand sides at the requested point. @qcode{"TolX"} specifies |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
42 ## the termination tolerance in the unknown variables, while |
17289
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17253
diff
changeset
|
43 ## @qcode{"TolFun"} is a tolerance for equations. Default is @code{1e-7} |
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17253
diff
changeset
|
44 ## for both @qcode{"TolX"} and @qcode{"TolFun"}. |
10201
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
45 ## |
17289
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17253
diff
changeset
|
46 ## If @qcode{"AutoScaling"} is on, the variables will be automatically scaled |
10711
fbd7843974fa
Periodic grammar check of documentation files to ensure common format.
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
47 ## according to the column norms of the (estimated) Jacobian. As a result, |
fbd7843974fa
Periodic grammar check of documentation files to ensure common format.
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
48 ## TolF becomes scaling-independent. By default, this option is off, because |
10201
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
49 ## it may sometimes deliver unexpected (though mathematically correct) results. |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
50 ## |
17289
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17253
diff
changeset
|
51 ## If @qcode{"Updating"} is @qcode{"on"}, the function will attempt to use |
17298
17be601bc783
doc: Periodic spellchecking of documentation.
Rik <rik@octave.org>
parents:
17289
diff
changeset
|
52 ## @nospell{Broyden} updates to update the Jacobian, in order to reduce the |
17be601bc783
doc: Periodic spellchecking of documentation.
Rik <rik@octave.org>
parents:
17289
diff
changeset
|
53 ## amount of Jacobian calculations. If your user function always calculates the |
17289
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17253
diff
changeset
|
54 ## Jacobian (regardless of number of output arguments), this option provides |
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17253
diff
changeset
|
55 ## no advantage and should be set to false. |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
56 ## |
17289
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17253
diff
changeset
|
57 ## @qcode{"ComplexEqn"} is @qcode{"on"}, @code{fsolve} will attempt to solve |
9065
8207b833557f
Cleanup documentation for arith.texi, linalg.texi, nonlin.texi
Rik <rdrider0-list@yahoo.com>
parents:
9051
diff
changeset
|
58 ## complex equations in complex variables, assuming that the equations possess a |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
59 ## complex derivative (i.e., are holomorphic). If this is not what you want, |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
60 ## 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
|
61 ## system. |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
62 ## |
8466 | 63 ## For description of the other options, see @code{optimset}. |
64 ## | |
65 ## On return, @var{fval} contains the value of the function @var{fcn} | |
66 ## evaluated at @var{x}, and @var{info} may be one of the following values: | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
67 ## |
8466 | 68 ## @table @asis |
69 ## @item 1 | |
10793
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
10711
diff
changeset
|
70 ## Converged to a solution point. Relative residual error is less than |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
71 ## specified by TolFun. |
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
72 ## |
8466 | 73 ## @item 2 |
74 ## Last relative step size was less that TolX. | |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
75 ## |
8466 | 76 ## @item 3 |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
77 ## Last relative decrease in residual was less than TolF. |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
78 ## |
8466 | 79 ## @item 0 |
80 ## Iteration limit exceeded. | |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
81 ## |
8466 | 82 ## @item -3 |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
83 ## The trust region radius became excessively small. |
8466 | 84 ## @end table |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
85 ## |
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
86 ## Note: If you only have a single nonlinear equation of one variable, using |
8466 | 87 ## @code{fzero} is usually a much better idea. |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
88 ## |
10793
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
10711
diff
changeset
|
89 ## Note about user-supplied Jacobians: |
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
10711
diff
changeset
|
90 ## As an inherent property of the algorithm, Jacobian is always requested for a |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
91 ## solution vector whose residual vector is already known, and it is the last |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8997
diff
changeset
|
92 ## accepted successful step. Often this will be one of the last two calls, but |
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8997
diff
changeset
|
93 ## not always. If the savings by reusing intermediate results from residual |
10793
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
10711
diff
changeset
|
94 ## calculation in Jacobian calculation are significant, the best strategy is to |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
95 ## employ OutputFcn: After a vector is evaluated for residuals, if OutputFcn is |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
96 ## called with that vector, then the intermediate results should be saved for |
10793
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
10711
diff
changeset
|
97 ## future Jacobian evaluation, and should be kept until a Jacobian evaluation |
16013
4efc0c1537df
Fix typo in the documentation of function fsolve
Rafael Laboissiere <rafael@laboissiere.net>
parents:
14868
diff
changeset
|
98 ## is requested or until OutputFcn is called with a different vector, in which |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8997
diff
changeset
|
99 ## case they should be dropped in favor of this most recent vector. A short |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
100 ## example how this can be achieved follows: |
8988
315828058e0d
minor doc fixes in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8986
diff
changeset
|
101 ## |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
102 ## @example |
9153
5247e89688e1
Eliminate most overfull errors when running texi2pdf for generating pdf documentation
Rik <rdrider0-list@yahoo.com>
parents:
9143
diff
changeset
|
103 ## function [fvec, fjac] = user_func (x, optimvalues, state) |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
104 ## persistent sav = [], sav0 = []; |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
105 ## if (nargin == 1) |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
106 ## ## evaluation call |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
107 ## if (nargout == 1) |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
108 ## sav0.x = x; # mark saved vector |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
109 ## ## calculate fvec, save results to sav0. |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
110 ## elseif (nargout == 2) |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
111 ## ## calculate fjac using sav. |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
112 ## endif |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
113 ## else |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
114 ## ## outputfcn call. |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
115 ## if (all (x == sav0.x)) |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
116 ## sav = sav0; |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
117 ## endif |
9209
923c7cb7f13f
Simplify TeXinfo files by eliminating redundant @iftex followed by @tex construction.
Rik <rdrider0-list@yahoo.com>
parents:
9207
diff
changeset
|
118 ## ## maybe output iteration status, etc. |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
119 ## endif |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
120 ## endfunction |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
121 ## |
14143
04dcbb8fb880
fsolve.m: Move @seealso to bottom in docstring to silence warning.
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
14138
diff
changeset
|
122 ## ## @dots{} |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
123 ## |
9153
5247e89688e1
Eliminate most overfull errors when running texi2pdf for generating pdf documentation
Rik <rdrider0-list@yahoo.com>
parents:
9143
diff
changeset
|
124 ## fsolve (@@user_func, x0, optimset ("OutputFcn", @@user_func, @dots{})) |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
125 ## @end example |
14143
04dcbb8fb880
fsolve.m: Move @seealso to bottom in docstring to silence warning.
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
14138
diff
changeset
|
126 ## @seealso{fzero, optimset} |
8466 | 127 ## @end deftypefn |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
128 |
13027
b9a89ca0fb75
prevent optimization functions from setting ans in workspace at startup
John W. Eaton <jwe@octave.org>
parents:
11587
diff
changeset
|
129 ## PKG_ADD: ## Discard result to avoid polluting workspace with ans at startup. |
b9a89ca0fb75
prevent optimization functions from setting ans in workspace at startup
John W. Eaton <jwe@octave.org>
parents:
11587
diff
changeset
|
130 ## PKG_ADD: [~] = __all_opts__ ("fsolve"); |
8648
ff61b53eb294
optimization: use PKG_ADD: comments instead of PKG_ADD file
John W. Eaton <jwe@octave.org>
parents:
8647
diff
changeset
|
131 |
8596
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
132 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
|
133 |
8596
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
134 ## Get default options if requested. |
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
135 if (nargin == 1 && ischar (fcn) && strcmp (fcn, 'defaults')) |
17253
7babcdb9bc13
Use ... instead of \ for line continuation marker.
Stefan Mahr <dac922@gmx.de>
parents:
16013
diff
changeset
|
136 x = optimset ("MaxIter", 400, "MaxFunEvals", Inf, ... |
10201
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
137 "Jacobian", "off", "TolX", 1e-7, "TolFun", 1e-7, |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
138 "OutputFcn", [], "Updating", "on", "FunValCheck", "off", |
10201
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
139 "ComplexEqn", "off", "FinDiffType", "central", |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
140 "TypicalX", [], "AutoScaling", "off"); |
8596
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
141 return; |
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
142 endif |
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
143 |
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
144 if (nargin < 2 || nargin > 3 || ! ismatrix (x0)) |
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
145 print_usage (); |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
146 endif |
8596
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
147 |
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
148 if (ischar (fcn)) |
9464
e598248a060d
safer str2func use in optim functions
Jaroslav Hajek <highegg@gmail.com>
parents:
9212
diff
changeset
|
149 fcn = str2func (fcn, "global"); |
10050
dc88a0b6472c
support old style jacobian for fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
150 elseif (iscell (fcn)) |
dc88a0b6472c
support old style jacobian for fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
151 fcn = @(x) make_fcn_jac (x, fcn{1}, fcn{2}); |
8306
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 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
154 xsiz = size (x0); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
155 n = numel (x0); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
156 |
8467
77b8d4aa2743
fsolve.m, fzero.m: style fixes; use strcmpi to compare options
John W. Eaton <jwe@octave.org>
parents:
8466
diff
changeset
|
157 has_jac = strcmpi (optimget (options, "Jacobian", "off"), "on"); |
9212
6feb27c38da1
support central differences in fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9209
diff
changeset
|
158 cdif = strcmpi (optimget (options, "FinDiffType", "central"), "central"); |
8596
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
159 maxiter = optimget (options, "MaxIter", 400); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
160 maxfev = optimget (options, "MaxFunEvals", Inf); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
161 outfcn = optimget (options, "OutputFcn"); |
8596
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
162 updating = strcmpi (optimget (options, "Updating", "on"), "on"); |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
163 complexeqn = strcmpi (optimget (options, "ComplexEqn", "off"), "on"); |
8590 | 164 |
10201
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
165 ## Get scaling matrix using the TypicalX option. If set to "auto", the |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
166 ## scaling matrix is estimated using the Jacobian. |
10201
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
167 typicalx = optimget (options, "TypicalX"); |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
168 if (isempty (typicalx)) |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
169 typicalx = ones (n, 1); |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
170 endif |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
171 autoscale = strcmpi (optimget (options, "AutoScaling", "off"), "on"); |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
172 if (! autoscale) |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
173 dg = 1 ./ typicalx; |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
174 endif |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
175 |
8467
77b8d4aa2743
fsolve.m, fzero.m: style fixes; use strcmpi to compare options
John W. Eaton <jwe@octave.org>
parents:
8466
diff
changeset
|
176 funvalchk = strcmpi (optimget (options, "FunValCheck", "off"), "on"); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
177 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
178 if (funvalchk) |
8604 | 179 ## Replace fcn with a guarded version. |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
180 fcn = @(x) guarded_eval (fcn, x, complexeqn); |
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 |
8466 | 183 ## These defaults are rather stringent. I think that normally, user |
184 ## prefers accuracy to performance. | |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
185 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
186 macheps = eps (class (x0)); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
187 |
10201
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
188 tolx = optimget (options, "TolX", 1e-7); |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
189 tolf = optimget (options, "TolFun", 1e-7); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
190 |
9623
bc0739d02724
update initial TR step for fsolve and fminunc
Jaroslav Hajek <highegg@gmail.com>
parents:
9464
diff
changeset
|
191 factor = 1; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
192 |
8467
77b8d4aa2743
fsolve.m, fzero.m: style fixes; use strcmpi to compare options
John W. Eaton <jwe@octave.org>
parents:
8466
diff
changeset
|
193 niter = 1; |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
194 nfev = 1; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
195 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
196 x = x0(:); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
197 info = 0; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
198 |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
199 ## Initial evaluation. |
8990
349d75161672
more cosmetic adjustments to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8988
diff
changeset
|
200 ## Handle arbitrary shapes of x and f and remember them. |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
201 fvec = fcn (reshape (x, xsiz)); |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
202 fsiz = size (fvec); |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
203 fvec = fvec(:); |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
204 fn = norm (fvec); |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
205 m = length (fvec); |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
206 n = length (x); |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
207 |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
208 if (! isempty (outfcn)) |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
209 optimvalues.iter = niter; |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
210 optimvalues.funccount = nfev; |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
211 optimvalues.fval = fn; |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
212 optimvalues.searchdirection = zeros (n, 1); |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
213 state = 'init'; |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
214 stop = outfcn (x, optimvalues, state); |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
215 if (stop) |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
216 info = -1; |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
217 break; |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
218 endif |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
219 endif |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
220 |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
221 nsuciter = 0; |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
222 |
8507 | 223 ## Outer loop. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
224 while (niter < maxiter && nfev < maxfev && ! info) |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
225 |
8507 | 226 ## Calculate function value and Jacobian (possibly via FD). |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
227 if (has_jac) |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
228 [fvec, fjac] = fcn (reshape (x, xsiz)); |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
229 ## If the Jacobian is sparse, disable Broyden updating. |
8592
dacfd030633a
handle sparse jacobians in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8590
diff
changeset
|
230 if (issparse (fjac)) |
dacfd030633a
handle sparse jacobians in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8590
diff
changeset
|
231 updating = false; |
dacfd030633a
handle sparse jacobians in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8590
diff
changeset
|
232 endif |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
233 fvec = fvec(:); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
234 nfev ++; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
235 else |
10201
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
236 fjac = __fdjac__ (fcn, reshape (x, xsiz), fvec, typicalx, cdif); |
9212
6feb27c38da1
support central differences in fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9209
diff
changeset
|
237 nfev += (1 + cdif) * length (x); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
238 endif |
8590 | 239 |
8763
5ce12bca4c51
update comments in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8693
diff
changeset
|
240 ## For square and overdetermined systems, we update a QR |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
241 ## factorization of the Jacobian to avoid solving a full system in each |
8590 | 242 ## step. In this case, we pass a triangular matrix to __dogleg__. |
243 useqr = updating && m >= n && n > 10; | |
244 | |
245 if (useqr) | |
8616
3d75d717cbe0
do not pivot by default in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8615
diff
changeset
|
246 ## 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
|
247 ## 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
|
248 ## 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
|
249 ## 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
|
250 ## when qr updating was used. |
3d75d717cbe0
do not pivot by default in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8615
diff
changeset
|
251 [q, r] = qr (fjac, 0); |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
252 endif |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
253 |
10201
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
254 if (autoscale) |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
255 ## Get column norms, use them as scaling factors. |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
256 jcn = norm (fjac, 'columns').'; |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
257 if (niter == 1) |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
258 dg = jcn; |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
259 dg(dg == 0) = 1; |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
260 else |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
261 ## Rescale adaptively. |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
262 ## FIXME: the original minpack used the following rescaling strategy: |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
263 ## dg = max (dg, jcn); |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
264 ## but it seems not good if we start with a bad guess yielding Jacobian |
10201
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
265 ## columns with large norms that later decrease, because the corresponding |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
266 ## variable will still be overscaled. So instead, we only give the old |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
267 ## scaling a small momentum, but do not honor it. |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
268 |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
269 dg = max (0.1*dg, jcn); |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
270 endif |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
271 endif |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
272 |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
273 if (niter == 1) |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
274 xn = norm (dg .* x); |
8590 | 275 ## FIXME: something better? |
9628
73e6ad869f08
further correct initial TR step strategy
Jaroslav Hajek <highegg@gmail.com>
parents:
9623
diff
changeset
|
276 delta = factor * max (xn, 1); |
8306
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 |
8990
349d75161672
more cosmetic adjustments to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8988
diff
changeset
|
279 ## It also seems that in the case of fast (and inhomogeneously) changing |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
280 ## Jacobian, the Broyden updates are of little use, so maybe we could |
8990
349d75161672
more cosmetic adjustments to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8988
diff
changeset
|
281 ## skip them if a big disproportional change is expected. The question is, |
349d75161672
more cosmetic adjustments to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8988
diff
changeset
|
282 ## of course, how to define the above terms :) |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
283 |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
284 lastratio = 0; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
285 nfail = 0; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
286 nsuc = 0; |
9207
25f50d2d76b3
improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9153
diff
changeset
|
287 decfac = 0.5; |
25f50d2d76b3
improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9153
diff
changeset
|
288 |
8507 | 289 ## Inner loop. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
290 while (niter <= maxiter && nfev < maxfev && ! info) |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
291 |
8590 | 292 ## Get trust-region model (dogleg) minimizer. |
293 if (useqr) | |
294 qtf = q'*fvec; | |
295 s = - __dogleg__ (r, qtf, dg, delta); | |
296 w = qtf + r * s; | |
297 else | |
298 s = - __dogleg__ (fjac, fvec, dg, delta); | |
299 w = fvec + fjac * s; | |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
300 endif |
8590 | 301 |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
302 sn = norm (dg .* s); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
303 if (niter == 1) |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
304 delta = min (delta, sn); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
305 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
306 |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
307 fvec1 = fcn (reshape (x + s, xsiz)) (:); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
308 fn1 = norm (fvec1); |
8590 | 309 nfev ++; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
310 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
311 if (fn1 < fn) |
8507 | 312 ## Scaled actual reduction. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
313 actred = 1 - (fn1/fn)^2; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
314 else |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
315 actred = -1; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
316 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
317 |
8507 | 318 ## Scaled predicted reduction, and ratio. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
319 t = norm (w); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
320 if (t < fn) |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
321 prered = 1 - (t/fn)^2; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
322 ratio = actred / prered; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
323 else |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
324 prered = 0; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
325 ratio = 0; |
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 ## Update delta. |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14386
diff
changeset
|
329 if (ratio < min (max (0.1, 0.8*lastratio), 0.9)) |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
330 nsuc = 0; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
331 nfail ++; |
9207
25f50d2d76b3
improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9153
diff
changeset
|
332 delta *= decfac; |
25f50d2d76b3
improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9153
diff
changeset
|
333 decfac ^= 1.4142; |
8590 | 334 if (delta <= 1e1*macheps*xn) |
8507 | 335 ## Trust region became uselessly small. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
336 info = -3; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
337 break; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
338 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
339 else |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
340 lastratio = ratio; |
9207
25f50d2d76b3
improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9153
diff
changeset
|
341 decfac = 0.5; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
342 nfail = 0; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
343 nsuc ++; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
344 if (abs (1-ratio) <= 0.1) |
9207
25f50d2d76b3
improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9153
diff
changeset
|
345 delta = 1.4142*sn; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
346 elseif (ratio >= 0.5 || nsuc > 1) |
9207
25f50d2d76b3
improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9153
diff
changeset
|
347 delta = max (delta, 1.4142*sn); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
348 endif |
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 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
351 if (ratio >= 1e-4) |
8507 | 352 ## Successful iteration. |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
353 x += s; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
354 xn = norm (dg .* x); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
355 fvec = fvec1; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
356 fn = fn1; |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
357 nsuciter ++; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
358 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
359 |
8857 | 360 niter ++; |
361 | |
8590 | 362 ## FIXME: should outputfcn be only called after a successful iteration? |
8615 | 363 if (! isempty (outfcn)) |
8590 | 364 optimvalues.iter = niter; |
365 optimvalues.funccount = nfev; | |
366 optimvalues.fval = fn; | |
367 optimvalues.searchdirection = s; | |
368 state = 'iter'; | |
369 stop = outfcn (x, optimvalues, state); | |
370 if (stop) | |
371 info = -1; | |
372 break; | |
373 endif | |
374 endif | |
375 | |
8466 | 376 ## Tests for termination conditions. A mysterious place, anything |
377 ## can happen if you change something here... | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
378 |
8466 | 379 ## The rule of thumb (which I'm not sure M*b is quite following) |
380 ## is that for a tolerance that depends on scaling, only 0 makes | |
381 ## sense as a default value. But 0 usually means uselessly long | |
382 ## iterations, so we need scaling-independent tolerances wherever | |
383 ## possible. | |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
384 |
8466 | 385 ## FIXME -- why tolf*n*xn? If abs (e) ~ abs(x) * eps is a vector |
386 ## of perturbations of x, then norm (fjac*e) <= eps*n*xn, i.e. by | |
387 ## 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
|
388 if (fn <= tolf*n*xn) |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
389 info = 1; |
10549 | 390 ## The following tests done only after successful step. |
9075 | 391 elseif (ratio >= 1e-4) |
8466 | 392 ## This one is classic. Note that we use scaled variables again, |
10549 | 393 ## but compare to scaled step, so nothing bad. |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
394 if (sn <= tolx*xn) |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
395 info = 2; |
8466 | 396 ## Again a classic one. It seems weird to use the same tolf |
10549 | 397 ## for two different tests, but that's what M*b manual appears |
398 ## to say. | |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
399 elseif (actred < tolf) |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
400 info = 3; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
401 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
402 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
403 |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
404 ## Criterion for recalculating Jacobian. |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
405 if (! updating || nfail == 2 || nsuciter < 2) |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
406 break; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
407 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
408 |
8507 | 409 ## Compute the scaled Broyden update. |
8590 | 410 if (useqr) |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
411 u = (fvec1 - q*w) / sn; |
8590 | 412 v = dg .* ((dg .* s) / sn); |
413 | |
414 ## Update the QR factorization. | |
415 [q, r] = qrupdate (q, r, u, v); | |
416 else | |
417 u = (fvec1 - w); | |
418 v = dg .* ((dg .* s) / sn); | |
419 | |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
420 ## update the Jacobian |
8590 | 421 fjac += u * v'; |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
422 endif |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
423 endwhile |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
424 endwhile |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
425 |
8507 | 426 ## Restore original shapes. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
427 x = reshape (x, xsiz); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
428 fvec = reshape (fvec, fsiz); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
429 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
430 output.iterations = niter; |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
431 output.successful = nsuciter; |
8590 | 432 output.funcCount = nfev; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
433 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
434 endfunction |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
435 |
8466 | 436 ## An assistant function that evaluates a function handle and checks for |
437 ## bad results. | |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
438 function [fx, jx] = guarded_eval (fun, x, complexeqn) |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
439 if (nargout > 1) |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
440 [fx, jx] = fun (x); |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
441 else |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
442 fx = fun (x); |
8997
187a9d9c2f04
simplifications to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8991
diff
changeset
|
443 jx = []; |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
444 endif |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
445 |
8997
187a9d9c2f04
simplifications to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8991
diff
changeset
|
446 if (! complexeqn && ! (isreal (fx) && isreal (jx))) |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
447 error ("fsolve:notreal", "fsolve: non-real value encountered"); |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14386
diff
changeset
|
448 elseif (complexeqn && ! (isnumeric (fx) && isnumeric (jx))) |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
449 error ("fsolve:notnum", "fsolve: non-numeric value encountered"); |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
450 elseif (any (isnan (fx(:)))) |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
451 error ("fsolve:isnan", "fsolve: NaN value encountered"); |
14386
59aab666f2bf
Extend "FunValCheck" option to optimization routines to detect Inf values.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
452 elseif (any (isinf (fx(:)))) |
59aab666f2bf
Extend "FunValCheck" option to optimization routines to detect Inf values.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
453 error ("fsolve:isinf", "fsolve: Inf value encountered"); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
454 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
455 endfunction |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
456 |
10050
dc88a0b6472c
support old style jacobian for fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
457 function [fx, jx] = make_fcn_jac (x, fcn, fjac) |
dc88a0b6472c
support old style jacobian for fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
458 fx = fcn (x); |
dc88a0b6472c
support old style jacobian for fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
459 if (nargout == 2) |
dc88a0b6472c
support old style jacobian for fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
460 jx = fjac (x); |
dc88a0b6472c
support old style jacobian for fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
461 endif |
dc88a0b6472c
support old style jacobian for fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
462 endfunction |
dc88a0b6472c
support old style jacobian for fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
463 |
17346
1c89599167a6
maint: End m-files with 1 blank line.
Rik <rik@octave.org>
parents:
17298
diff
changeset
|
464 |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
465 %!function retval = __f (p) |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
466 %! x = p(1); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
467 %! y = p(2); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
468 %! z = p(3); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
469 %! retval = zeros (3, 1); |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
470 %! retval(1) = sin (x) + y^2 + log (z) - 7; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
471 %! retval(2) = 3*x + 2^y -z^3 + 1; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
472 %! retval(3) = x + y + z - 5; |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
473 %!endfunction |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
474 %!test |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
475 %! x_opt = [ 0.599054; |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
476 %! 2.395931; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
477 %! 2.005014 ]; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
478 %! tol = 1.0e-5; |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
479 %! [x, fval, info] = fsolve (@__f, [ 0.5; 2.0; 2.5 ]); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
480 %! assert (info > 0); |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
481 %! assert (norm (x - x_opt, Inf) < tol); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
482 %! assert (norm (fval) < tol); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
483 |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
484 %!function retval = __f (p) |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
485 %! x = p(1); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
486 %! y = p(2); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
487 %! z = p(3); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
488 %! w = p(4); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
489 %! retval = zeros (4, 1); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
490 %! retval(1) = 3*x + 4*y + exp (z + w) - 1.007; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
491 %! retval(2) = 6*x - 4*y + exp (3*z + w) - 11; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
492 %! 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
|
493 %! retval(4) = x^2 + 2*y^3 + z - w - 4; |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
494 %!endfunction |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
495 %!test |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
496 %! x_opt = [ -0.767297326653401, 0.590671081117440, ... |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
497 %! 1.47190018629642, -1.52719341133957 ]; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
498 %! tol = 1.0e-5; |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
499 %! [x, fval, info] = fsolve (@__f, [-1, 1, 2, -1]); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
500 %! assert (info > 0); |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
501 %! assert (norm (x - x_opt, Inf) < tol); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
502 %! assert (norm (fval) < tol); |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
503 |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
504 %!function retval = __f (p) |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
505 %! x = p(1); |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
506 %! y = p(2); |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
507 %! z = p(3); |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
508 %! retval = zeros (3, 1); |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
509 %! retval(1) = sin (x) + y^2 + log (z) - 7; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
510 %! retval(2) = 3*x + 2^y -z^3 + 1; |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
511 %! retval(3) = x + y + z - 5; |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
512 %! retval(4) = x*x + y - z*log (z) - 1.36; |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
513 %!endfunction |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
514 %!test |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
515 %! x_opt = [ 0.599054; |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
516 %! 2.395931; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
517 %! 2.005014 ]; |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
518 %! tol = 1.0e-5; |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
519 %! [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
|
520 %! assert (info > 0); |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
521 %! assert (norm (x - x_opt, Inf) < tol); |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
522 %! assert (norm (fval) < tol); |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
523 |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
524 %!function retval = __f (p) |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
525 %! x = p(1); |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
526 %! y = p(2); |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
527 %! z = p(3); |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
528 %! retval = zeros (3, 1); |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
529 %! retval(1) = sin (x) + y^2 + log (z) - 7; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
530 %! retval(2) = 3*x + 2^y -z^3 + 1; |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
531 %! retval(3) = x + y + z - 5; |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
532 %!endfunction |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
533 %!test |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
534 %! x_opt = [ 0.599054; |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
535 %! 2.395931; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
536 %! 2.005014 ]; |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
537 %! tol = 1.0e-5; |
8590 | 538 %! opt = optimset ("Updating", "qrp"); |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
539 %! [x, fval, info] = fsolve (@__f, [ 0.5; 2.0; 2.5 ], opt); |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
540 %! assert (info > 0); |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
541 %! assert (norm (x - x_opt, Inf) < tol); |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
542 %! assert (norm (fval) < tol); |
8590 | 543 |
544 %!test | |
545 %! b0 = 3; | |
546 %! a0 = 0.2; | |
547 %! x = 0:.5:5; | |
548 %! noise = 1e-5 * sin (100*x); | |
549 %! y = exp (-a0*x) + b0 + noise; | |
550 %! c_opt = [a0, b0]; | |
551 %! tol = 1e-5; | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
552 %! |
8590 | 553 %! [c, fval, info, output] = fsolve (@(c) (exp(-c(1)*x) + c(2) - y), [0, 0]); |
554 %! assert (info > 0); | |
555 %! assert (norm (c - c_opt, Inf) < tol); | |
556 %! assert (norm (fval) < norm (noise)); | |
557 | |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
558 %!function y = cfun (x) |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
559 %! 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
|
560 %! 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
|
561 %! y(3) = x(1) * x(2) - x(3)^2 + (3+2i); |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
562 %!endfunction |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
563 |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
564 %!test |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
565 %! x_opt = [-1+i, 1-i, 2+i]; |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
566 %! x = [i, 1, 1+i]; |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
567 %! |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
568 %! [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
|
569 %! tol = 1e-5; |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
570 %! assert (norm (f) < tol); |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
571 %! assert (norm (x - x_opt, Inf) < tol); |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
572 |
9899
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
573 ## Solve the double dogleg trust-region least-squares problem: |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
574 ## Minimize norm(r*x-b) subject to the constraint norm(d.*x) <= delta, |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
575 ## x being a convex combination of the gauss-newton and scaled gradient. |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
576 |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
577 ## TODO: error checks |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
578 ## TODO: handle singularity, or leave it up to mldivide? |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
579 |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
580 function x = __dogleg__ (r, b, d, delta) |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
581 ## Get Gauss-Newton direction. |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
582 x = r \ b; |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
583 xn = norm (d .* x); |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
584 if (xn > delta) |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
585 ## GN is too big, get scaled gradient. |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
586 s = (r' * b) ./ d; |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
587 sn = norm (s); |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
588 if (sn > 0) |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
589 ## Normalize and rescale. |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
590 s = (s / sn) ./ d; |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
591 ## Get the line minimizer in s direction. |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
592 tn = norm (r*s); |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
593 snm = (sn / tn) / tn; |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
594 if (snm < delta) |
10549 | 595 ## Get the dogleg path minimizer. |
9899
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
596 bn = norm (b); |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
597 dxn = delta/xn; snmd = snm/delta; |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
598 t = (bn/sn) * (bn/xn) * snmd; |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
599 t -= dxn * snmd^2 - sqrt ((t-dxn)^2 + (1-dxn^2)*(1-snmd^2)); |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
600 alpha = dxn*(1-snmd^2) / t; |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
601 else |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
602 alpha = 0; |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
603 endif |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
604 else |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
605 alpha = delta / xn; |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
606 snm = 0; |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
607 endif |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
608 ## Form the appropriate convex combination. |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
609 x = alpha * x + ((1-alpha) * min (snm, delta)) * s; |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
610 endif |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
611 endfunction |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9724
diff
changeset
|
612 |