annotate scripts/optimization/fminunc.m @ 20695:6faaab833605

fminunc.m: Clean up function to meet Octave coding standards. * fminunc.m: Rename xsiz variable to xsz. Use double quotes for strings. Use in-place accumulation operator +=1 instead of ++ for performance. Use '__' prefix and postfix for internal functions. Use two spaces between end of sentence period and start of next sentence. Only calculate third output argument if requested.
author Rik <rik@octave.org>
date Wed, 11 Nov 2015 16:58:18 -0800
parents a7dbc4fc3732
children 5b7643257978
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
19697
4197fc428c7d maint: Update copyright notices for 2015.
John W. Eaton <jwe@octave.org>
parents: 18816
diff changeset
1 ## Copyright (C) 2008-2015 VZLU Prague, a.s.
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
2 ##
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
3 ## This file is part of Octave.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
4 ##
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
5 ## Octave is free software; you can redistribute it and/or modify it
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
6 ## under the terms of the GNU General Public License as published by
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
7 ## the Free Software Foundation; either version 3 of the License, or (at
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
8 ## your option) any later version.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
9 ##
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
10 ## Octave is distributed in the hope that it will be useful, but
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
11 ## WITHOUT ANY WARRANTY; without even the implied warranty of
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
12 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
13 ## General Public License for more details.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
14 ##
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
15 ## You should have received a copy of the GNU General Public License
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
16 ## along with Octave; see the file COPYING. If not, see
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
17 ## <http://www.gnu.org/licenses/>.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
18 ##
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
19 ## Author: Jaroslav Hajek <highegg@gmail.com>
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
20
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
21 ## -*- texinfo -*-
10687
a8ce6bdecce5 Improve documentation strings.
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
22 ## @deftypefn {Function File} {} fminunc (@var{fcn}, @var{x0})
a8ce6bdecce5 Improve documentation strings.
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
23 ## @deftypefnx {Function File} {} fminunc (@var{fcn}, @var{x0}, @var{options})
18608
cdcf66f4e244 fminunc.m: Fix typo in header documentation (bug #42011).
Massimiliano Fasi <mogrob.sanit@gmail.com>
parents: 17744
diff changeset
24 ## @deftypefnx {Function File} {[@var{x}, @var{fval}, @var{info}, @var{output}, @var{grad}, @var{hess}] =} fminunc (@var{fcn}, @dots{})
10793
be55736a0783 Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents: 10687
diff changeset
25 ## Solve an unconstrained optimization problem defined by the function
be55736a0783 Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents: 10687
diff changeset
26 ## @var{fcn}.
14895
e0525ecf156e Add new function fminsearch.m
Andy Adler <andy@analyti.ca>
parents: 14868
diff changeset
27 ##
20165
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
28 ## @var{fcn} should accept a vector (array) defining the unknown variables, and
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
29 ## return the objective function value, optionally with gradient.
18609
923614060f1d fminunc.m: Improve documentation.
Rik <rik@octave.org>
parents: 18608
diff changeset
30 ## @code{fminunc} attempts to determine a vector @var{x} such that
20165
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
31 ## @code{@var{fcn} (@var{x})} is a local minimum.
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
32 ##
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
33 ## @var{x0} determines a starting guess. The shape of @var{x0} is preserved in
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
34 ## all calls to @var{fcn}, but otherwise is treated as a column vector.
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
35 ##
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
36 ## @var{options} is a structure specifying additional options. Currently,
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
37 ## @code{fminunc} recognizes these options:
17281
bc924baa2c4e doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents: 14895
diff changeset
38 ## @qcode{"FunValCheck"}, @qcode{"OutputFcn"}, @qcode{"TolX"},
bc924baa2c4e doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents: 14895
diff changeset
39 ## @qcode{"TolFun"}, @qcode{"MaxIter"}, @qcode{"MaxFunEvals"},
20165
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
40 ## @qcode{"GradObj"}, @qcode{"FinDiffType"}, @qcode{"TypicalX"},
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
41 ## @qcode{"AutoScaling"}.
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
42 ##
20165
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
43 ## If @qcode{"GradObj"} is @qcode{"on"}, it specifies that @var{fcn}, when
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
44 ## called with two output arguments, also returns the Jacobian matrix of partial
20165
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
45 ## first derivatives at the requested point. @code{TolX} specifies the
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
46 ## termination tolerance for the unknown variables @var{x}, while @code{TolFun}
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
47 ## is a tolerance for the objective function value @var{fval}. The default is
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
48 ## @code{1e-7} for both options.
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
49 ##
18609
923614060f1d fminunc.m: Improve documentation.
Rik <rik@octave.org>
parents: 18608
diff changeset
50 ## For a description of the other options, see @code{optimset}.
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
51 ##
18609
923614060f1d fminunc.m: Improve documentation.
Rik <rik@octave.org>
parents: 18608
diff changeset
52 ## On return, @var{x} is the location of the minimum and @var{fval} contains
20165
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
53 ## the value of the objective function at @var{x}.
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
54 ##
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
55 ## @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
56 ##
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
57 ## @table @asis
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
58 ## @item 1
10793
be55736a0783 Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents: 10687
diff changeset
59 ## Converged to a solution point. Relative gradient error is less than
18609
923614060f1d fminunc.m: Improve documentation.
Rik <rik@octave.org>
parents: 18608
diff changeset
60 ## specified by @code{TolFun}.
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10793
diff changeset
61 ##
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
62 ## @item 2
18609
923614060f1d fminunc.m: Improve documentation.
Rik <rik@octave.org>
parents: 18608
diff changeset
63 ## Last relative step size was less than @code{TolX}.
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10793
diff changeset
64 ##
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
65 ## @item 3
18609
923614060f1d fminunc.m: Improve documentation.
Rik <rik@octave.org>
parents: 18608
diff changeset
66 ## Last relative change in function value was less than @code{TolFun}.
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10793
diff changeset
67 ##
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
68 ## @item 0
18816
e275d15c27b5 doc: Periodic spellcheck of documentation.
Rik <rik@octave.org>
parents: 18609
diff changeset
69 ## Iteration limit exceeded---either maximum number of algorithm iterations
18609
923614060f1d fminunc.m: Improve documentation.
Rik <rik@octave.org>
parents: 18608
diff changeset
70 ## @code{MaxIter} or maximum number of function evaluations @code{MaxFunEvals}.
923614060f1d fminunc.m: Improve documentation.
Rik <rik@octave.org>
parents: 18608
diff changeset
71 ##
923614060f1d fminunc.m: Improve documentation.
Rik <rik@octave.org>
parents: 18608
diff changeset
72 ## @item -1
18816
e275d15c27b5 doc: Periodic spellcheck of documentation.
Rik <rik@octave.org>
parents: 18609
diff changeset
73 ## Algorithm terminated by @code{OutputFcn}.
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10793
diff changeset
74 ##
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
75 ## @item -3
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
76 ## The trust region radius became excessively small.
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
77 ## @end table
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
78 ##
18609
923614060f1d fminunc.m: Improve documentation.
Rik <rik@octave.org>
parents: 18608
diff changeset
79 ## Optionally, @code{fminunc} can return a structure with convergence statistics
923614060f1d fminunc.m: Improve documentation.
Rik <rik@octave.org>
parents: 18608
diff changeset
80 ## (@var{output}), the output gradient (@var{grad}) at the solution @var{x},
923614060f1d fminunc.m: Improve documentation.
Rik <rik@octave.org>
parents: 18608
diff changeset
81 ## and approximate Hessian (@var{hess}) at the solution @var{x}.
9627
5bcfa0b346e8 fix extra outputs in fminunc
Jaroslav Hajek <highegg@gmail.com>
parents: 9623
diff changeset
82 ##
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
83 ## Application Notes: If the objective function is a single nonlinear equation
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
84 ## of one variable then using @code{fminbnd} is usually a better choice.
20165
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
85 ##
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
86 ## The algorithm used by @code{fminsearch} is a gradient search which depends
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
87 ## on the objective function being differentiable. If the function has
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
88 ## discontinuities it may be better to use a derivative-free algorithm such as
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
89 ## @code{fminsearch}.
14895
e0525ecf156e Add new function fminsearch.m
Andy Adler <andy@analyti.ca>
parents: 14868
diff changeset
90 ## @seealso{fminbnd, fminsearch, optimset}
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
91 ## @end deftypefn
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
92
13027
b9a89ca0fb75 prevent optimization functions from setting ans in workspace at startup
John W. Eaton <jwe@octave.org>
parents: 12578
diff changeset
93 ## 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: 12578
diff changeset
94 ## PKG_ADD: [~] = __all_opts__ ("fminunc");
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
95
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
96 function [x, fval, info, output, grad, hess] = fminunc (fcn, x0, options = struct ())
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
97
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
98 ## Get default options if requested.
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
99 if (nargin == 1 && strcmp (fcn, "defaults"))
14552
86854d032a37 maint: miscellaneous style fixes for .m files
John W. Eaton <jwe@octave.org>
parents: 14386
diff changeset
100 x = optimset ("MaxIter", 400, "MaxFunEvals", Inf,
86854d032a37 maint: miscellaneous style fixes for .m files
John W. Eaton <jwe@octave.org>
parents: 14386
diff changeset
101 "GradObj", "off", "TolX", 1e-7, "TolFun", 1e-7,
86854d032a37 maint: miscellaneous style fixes for .m files
John W. Eaton <jwe@octave.org>
parents: 14386
diff changeset
102 "OutputFcn", [], "FunValCheck", "off",
86854d032a37 maint: miscellaneous style fixes for .m files
John W. Eaton <jwe@octave.org>
parents: 14386
diff changeset
103 "FinDiffType", "central",
86854d032a37 maint: miscellaneous style fixes for .m files
John W. Eaton <jwe@octave.org>
parents: 14386
diff changeset
104 "TypicalX", [], "AutoScaling", "off");
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
105 return;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
106 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
107
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
108 if (nargin < 2 || nargin > 3 || ! ismatrix (x0))
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
109 print_usage ();
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
110 endif
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
111
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
112 if (ischar (fcn))
9464
e598248a060d safer str2func use in optim functions
Jaroslav Hajek <highegg@gmail.com>
parents: 9212
diff changeset
113 fcn = str2func (fcn, "global");
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
114 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
115
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
116 xsz = size (x0);
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
117 n = numel (x0);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
118
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
119 has_grad = strcmpi (optimget (options, "GradObj", "off"), "on");
9212
6feb27c38da1 support central differences in fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9207
diff changeset
120 cdif = strcmpi (optimget (options, "FinDiffType", "central"), "central");
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
121 maxiter = optimget (options, "MaxIter", 400);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
122 maxfev = optimget (options, "MaxFunEvals", Inf);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
123 outfcn = optimget (options, "OutputFcn");
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
124
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
125 ## Get scaling matrix using the TypicalX option. If set to "auto", the
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
126 ## 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: 9899
diff changeset
127 typicalx = optimget (options, "TypicalX");
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
128 if (isempty (typicalx))
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
129 typicalx = ones (n, 1);
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
130 endif
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
131 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: 9899
diff changeset
132 if (! autoscale)
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
133 dg = 1 ./ typicalx;
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
134 endif
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
135
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
136 funvalchk = strcmpi (optimget (options, "FunValCheck", "off"), "on");
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
137
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
138 if (funvalchk)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
139 ## Replace fcn with a guarded version.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
140 fcn = @(x) guarded_eval (fcn, x);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
141 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
142
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
143 ## These defaults are rather stringent. I think that normally, user
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
144 ## prefers accuracy to performance.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
145
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
146 macheps = eps (class (x0));
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
147
10201
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
148 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: 9899
diff changeset
149 tolf = optimget (options, "TolFun", 1e-7);
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
150
9623
bc0739d02724 update initial TR step for fsolve and fminunc
Jaroslav Hajek <highegg@gmail.com>
parents: 9464
diff changeset
151 factor = 0.1;
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
152 ## FIXME: TypicalX corresponds to user scaling (???)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
153 autodg = true;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
154
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
155 niter = 1;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
156 nfev = 0;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
157
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
158 x = x0(:);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
159 info = 0;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
160
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
161 ## Initial evaluation.
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
162 fval = fcn (reshape (x, xsz));
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
163 n = length (x);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
164
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
165 if (! isempty (outfcn))
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
166 optimvalues.iter = niter;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
167 optimvalues.funccount = nfev;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
168 optimvalues.fval = fval;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
169 optimvalues.searchdirection = zeros (n, 1);
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
170 state = "init";
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
171 stop = outfcn (x, optimvalues, state);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
172 if (stop)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
173 info = -1;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
174 break;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
175 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
176 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
177
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
178 nsuciter = 0;
9199
399884c9d4a1 import the step adaptation strategy from fsolve to fminunc
Jaroslav Hajek <highegg@gmail.com>
parents: 9084
diff changeset
179 lastratio = 0;
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
180
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
181 grad = [];
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
182
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
183 ## Outer loop.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
184 while (niter < maxiter && nfev < maxfev && ! info)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
185
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
186 grad0 = grad;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
187
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
188 ## Calculate function value and gradient (possibly via FD).
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
189 if (has_grad)
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
190 [fval, grad] = fcn (reshape (x, xsz));
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
191 grad = grad(:);
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
192 nfev += 1;
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
193 else
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
194 grad = __fdjac__ (fcn, reshape (x, xsz), fval, typicalx, cdif)(:);
9212
6feb27c38da1 support central differences in fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9207
diff changeset
195 nfev += (1 + cdif) * length (x);
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
196 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
197
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
198 if (niter == 1)
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
199 ## Initialize by identity matrix.
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
200 hesr = eye (n);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
201 else
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
202 ## Use the damped BFGS formula.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
203 y = grad - grad0;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
204 sBs = sumsq (w);
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
205 Bs = hesr' * w;
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
206 sy = y' * s;
9633
ecc2c556f844 simplify damped BFGS formula
Jaroslav Hajek <highegg@gmail.com>
parents: 9631
diff changeset
207 theta = 0.8 / max (1 - sy / sBs, 0.8);
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
208 r = theta * y + (1-theta) * Bs;
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
209 hesr = cholupdate (hesr, r / sqrt (s' * r), "+");
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
210 [hesr, info] = cholupdate (hesr, Bs / sqrt (sBs), "-");
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
211 if (info)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
212 hesr = eye (n);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
213 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
214 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
215
10201
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
216 if (autoscale)
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
217 ## Second derivatives approximate the Hessian.
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
218 d2f = norm (hesr, "columns").';
10201
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
219 if (niter == 1)
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
220 dg = d2f;
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
221 else
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
222 ## FIXME: maybe fixed lower and upper bounds?
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
223 dg = max (0.1*dg, d2f);
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
224 endif
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
225 endif
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
226
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
227 if (niter == 1)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
228 xn = norm (dg .* x);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
229 ## FIXME: something better?
9628
73e6ad869f08 further correct initial TR step strategy
Jaroslav Hajek <highegg@gmail.com>
parents: 9627
diff changeset
230 delta = factor * max (xn, 1);
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
231 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
232
18609
923614060f1d fminunc.m: Improve documentation.
Rik <rik@octave.org>
parents: 18608
diff changeset
233 ## FIXME: why tolf*n*xn? If abs (e) ~ abs(x) * eps is a vector
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
234 ## of perturbations of x, then norm (hesr*e) <= eps*xn, i.e. by
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
235 ## tolf ~ eps we demand as much accuracy as we can expect.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
236 if (norm (grad) <= tolf*n*xn)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
237 info = 1;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
238 break;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
239 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
240
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
241 suc = false;
9207
25f50d2d76b3 improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9199
diff changeset
242 decfac = 0.5;
25f50d2d76b3 improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9199
diff changeset
243
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
244 ## Inner loop.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
245 while (! suc && niter <= maxiter && nfev < maxfev && ! info)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
246
9631
00958d0c4e3c split __dogleg__ > __doglegm__
Jaroslav Hajek <highegg@gmail.com>
parents: 9628
diff changeset
247 s = - __doglegm__ (hesr, grad, dg, delta);
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
248
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
249 sn = norm (dg .* s);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
250 if (niter == 1)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
251 delta = min (delta, sn);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
252 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
253
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
254 fval1 = fcn (reshape (x + s, xsz))(:);
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
255 nfev += 1;
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
256
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
257 if (fval1 < fval)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
258 ## Scaled actual reduction.
19833
9fc020886ae9 maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19697
diff changeset
259 actred = (fval - fval1) / (abs (fval1) + abs (fval));
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
260 else
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
261 actred = -1;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
262 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
263
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
264 w = hesr * s;
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
265 ## Scaled predicted reduction, and ratio.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
266 t = 1/2 * sumsq (w) + grad'*s;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
267 if (t < 0)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
268 prered = -t/(abs (fval) + abs (fval + t));
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
269 ratio = actred / prered;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
270 else
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
271 prered = 0;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
272 ratio = 0;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
273 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
274
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
275 ## Update delta.
14868
5d3a684236b0 maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents: 14552
diff changeset
276 if (ratio < min (max (0.1, 0.8*lastratio), 0.9))
9207
25f50d2d76b3 improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9199
diff changeset
277 delta *= decfac;
25f50d2d76b3 improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9199
diff changeset
278 decfac ^= 1.4142;
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
279 if (delta <= 10*macheps*xn)
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
280 ## Trust region became uselessly small.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
281 info = -3;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
282 break;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
283 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
284 else
9199
399884c9d4a1 import the step adaptation strategy from fsolve to fminunc
Jaroslav Hajek <highegg@gmail.com>
parents: 9084
diff changeset
285 lastratio = ratio;
9207
25f50d2d76b3 improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9199
diff changeset
286 decfac = 0.5;
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
287 if (abs (1-ratio) <= 0.1)
9207
25f50d2d76b3 improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9199
diff changeset
288 delta = 1.4142*sn;
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
289 elseif (ratio >= 0.5)
9207
25f50d2d76b3 improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9199
diff changeset
290 delta = max (delta, 1.4142*sn);
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
291 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
292 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
293
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
294 if (ratio >= 1e-4)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
295 ## Successful iteration.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
296 x += s;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
297 xn = norm (dg .* x);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
298 fval = fval1;
18609
923614060f1d fminunc.m: Improve documentation.
Rik <rik@octave.org>
parents: 18608
diff changeset
299 nsuciter++;
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
300 suc = true;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
301 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
302
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
303 niter += 1;
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
304
18609
923614060f1d fminunc.m: Improve documentation.
Rik <rik@octave.org>
parents: 18608
diff changeset
305 ## FIXME: should outputfcn be called only after a successful iteration?
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
306 if (! isempty (outfcn))
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
307 optimvalues.iter = niter;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
308 optimvalues.funccount = nfev;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
309 optimvalues.fval = fval;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
310 optimvalues.searchdirection = s;
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
311 state = "iter";
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
312 stop = outfcn (x, optimvalues, state);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
313 if (stop)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
314 info = -1;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
315 break;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
316 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
317 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
318
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
319 ## Tests for termination conditions. A mysterious place, anything
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
320 ## 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
321
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
322 ## The rule of thumb (which I'm not sure M*b is quite following)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
323 ## is that for a tolerance that depends on scaling, only 0 makes
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
324 ## sense as a default value. But 0 usually means uselessly long
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
325 ## iterations, so we need scaling-independent tolerances wherever
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
326 ## possible.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
327
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
328 ## The following tests done only after successful step.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
329 if (ratio >= 1e-4)
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
330 ## This one is classic. Note that we use scaled variables again,
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
331 ## but compare to scaled step, so nothing bad.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
332 if (sn <= tolx*xn)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
333 info = 2;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
334 ## Again a classic one.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
335 elseif (actred < tolf)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
336 info = 3;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
337 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
338 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
339
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
340 endwhile
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
341 endwhile
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
342
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
343 ## When info != 1, recalculate the gradient and Hessian using the final x.
20694
a7dbc4fc3732 fminunc.m: Correctly calculate gradient and hessian even when fminunc fails (bug #40101)
José Vallet <jose.vallet@aalto.fi>
parents: 20165
diff changeset
344 if (nargout > 4 && (info == -1 || info == 2 || info == 3))
a7dbc4fc3732 fminunc.m: Correctly calculate gradient and hessian even when fminunc fails (bug #40101)
José Vallet <jose.vallet@aalto.fi>
parents: 20165
diff changeset
345 grad0 = grad;
a7dbc4fc3732 fminunc.m: Correctly calculate gradient and hessian even when fminunc fails (bug #40101)
José Vallet <jose.vallet@aalto.fi>
parents: 20165
diff changeset
346 if (has_grad)
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
347 [fval, grad] = fcn (reshape (x, xsz));
20694
a7dbc4fc3732 fminunc.m: Correctly calculate gradient and hessian even when fminunc fails (bug #40101)
José Vallet <jose.vallet@aalto.fi>
parents: 20165
diff changeset
348 grad = grad(:);
a7dbc4fc3732 fminunc.m: Correctly calculate gradient and hessian even when fminunc fails (bug #40101)
José Vallet <jose.vallet@aalto.fi>
parents: 20165
diff changeset
349 else
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
350 grad = __fdjac__ (fcn, reshape (x, xsz), fval, typicalx, cdif)(:);
20694
a7dbc4fc3732 fminunc.m: Correctly calculate gradient and hessian even when fminunc fails (bug #40101)
José Vallet <jose.vallet@aalto.fi>
parents: 20165
diff changeset
351 endif
a7dbc4fc3732 fminunc.m: Correctly calculate gradient and hessian even when fminunc fails (bug #40101)
José Vallet <jose.vallet@aalto.fi>
parents: 20165
diff changeset
352
a7dbc4fc3732 fminunc.m: Correctly calculate gradient and hessian even when fminunc fails (bug #40101)
José Vallet <jose.vallet@aalto.fi>
parents: 20165
diff changeset
353 if (nargout > 5)
a7dbc4fc3732 fminunc.m: Correctly calculate gradient and hessian even when fminunc fails (bug #40101)
José Vallet <jose.vallet@aalto.fi>
parents: 20165
diff changeset
354 ## Use the damped BFGS formula.
a7dbc4fc3732 fminunc.m: Correctly calculate gradient and hessian even when fminunc fails (bug #40101)
José Vallet <jose.vallet@aalto.fi>
parents: 20165
diff changeset
355 y = grad - grad0;
a7dbc4fc3732 fminunc.m: Correctly calculate gradient and hessian even when fminunc fails (bug #40101)
José Vallet <jose.vallet@aalto.fi>
parents: 20165
diff changeset
356 sBs = sumsq (w);
a7dbc4fc3732 fminunc.m: Correctly calculate gradient and hessian even when fminunc fails (bug #40101)
José Vallet <jose.vallet@aalto.fi>
parents: 20165
diff changeset
357 Bs = hesr' * w;
a7dbc4fc3732 fminunc.m: Correctly calculate gradient and hessian even when fminunc fails (bug #40101)
José Vallet <jose.vallet@aalto.fi>
parents: 20165
diff changeset
358 sy = y' * s;
a7dbc4fc3732 fminunc.m: Correctly calculate gradient and hessian even when fminunc fails (bug #40101)
José Vallet <jose.vallet@aalto.fi>
parents: 20165
diff changeset
359 theta = 0.8 / max (1 - sy / sBs, 0.8);
a7dbc4fc3732 fminunc.m: Correctly calculate gradient and hessian even when fminunc fails (bug #40101)
José Vallet <jose.vallet@aalto.fi>
parents: 20165
diff changeset
360 r = theta * y + (1-theta) * Bs;
a7dbc4fc3732 fminunc.m: Correctly calculate gradient and hessian even when fminunc fails (bug #40101)
José Vallet <jose.vallet@aalto.fi>
parents: 20165
diff changeset
361 hesr = cholupdate (hesr, r / sqrt (s' * r), "+");
a7dbc4fc3732 fminunc.m: Correctly calculate gradient and hessian even when fminunc fails (bug #40101)
José Vallet <jose.vallet@aalto.fi>
parents: 20165
diff changeset
362 hesr = cholupdate (hesr, Bs / sqrt (sBs), "-");
a7dbc4fc3732 fminunc.m: Correctly calculate gradient and hessian even when fminunc fails (bug #40101)
José Vallet <jose.vallet@aalto.fi>
parents: 20165
diff changeset
363 endif
a7dbc4fc3732 fminunc.m: Correctly calculate gradient and hessian even when fminunc fails (bug #40101)
José Vallet <jose.vallet@aalto.fi>
parents: 20165
diff changeset
364 ## Return the gradient in the same shape as x
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
365 grad = reshape (grad, xsz);
20694
a7dbc4fc3732 fminunc.m: Correctly calculate gradient and hessian even when fminunc fails (bug #40101)
José Vallet <jose.vallet@aalto.fi>
parents: 20165
diff changeset
366 endif
a7dbc4fc3732 fminunc.m: Correctly calculate gradient and hessian even when fminunc fails (bug #40101)
José Vallet <jose.vallet@aalto.fi>
parents: 20165
diff changeset
367
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
368 ## Restore original shapes.
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
369 x = reshape (x, xsz);
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
370
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
371 if (nargout > 3)
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
372 output.iterations = niter;
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
373 output.successful = nsuciter;
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
374 output.funcCount = nfev;
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
375 endif
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
376
9627
5bcfa0b346e8 fix extra outputs in fminunc
Jaroslav Hajek <highegg@gmail.com>
parents: 9623
diff changeset
377 if (nargout > 5)
5bcfa0b346e8 fix extra outputs in fminunc
Jaroslav Hajek <highegg@gmail.com>
parents: 9623
diff changeset
378 hess = hesr'*hesr;
5bcfa0b346e8 fix extra outputs in fminunc
Jaroslav Hajek <highegg@gmail.com>
parents: 9623
diff changeset
379 endif
5bcfa0b346e8 fix extra outputs in fminunc
Jaroslav Hajek <highegg@gmail.com>
parents: 9623
diff changeset
380
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
381 endfunction
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
382
18609
923614060f1d fminunc.m: Improve documentation.
Rik <rik@octave.org>
parents: 18608
diff changeset
383 ## A helper function that evaluates a function and checks for bad results.
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
384 function [fx, gx] = guarded_eval (fun, x)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
385 if (nargout > 1)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
386 [fx, gx] = fun (x);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
387 else
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
388 fx = fun (x);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
389 gx = [];
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
390 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
391
14383
07c55bceca23 Fix guarded_eval() subfunction in fminunc (bug #35534).
Olaf Till <olaf.till@uni-jena.de>
parents: 14138
diff changeset
392 if (! (isreal (fx) && isreal (gx)))
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
393 error ("fminunc:notreal", "fminunc: non-real value encountered");
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
394 elseif (any (isnan (fx(:))))
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
395 error ("fminunc:isnan", "fminunc: NaN value encountered");
14386
59aab666f2bf Extend "FunValCheck" option to optimization routines to detect Inf values.
Rik <octave@nomad.inbox5.com>
parents: 14383
diff changeset
396 elseif (any (isinf (fx(:))))
59aab666f2bf Extend "FunValCheck" option to optimization routines to detect Inf values.
Rik <octave@nomad.inbox5.com>
parents: 14383
diff changeset
397 error ("fminunc:isinf", "fminunc: Inf value encountered");
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
398 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
399 endfunction
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
400
14386
59aab666f2bf Extend "FunValCheck" option to optimization routines to detect Inf values.
Rik <octave@nomad.inbox5.com>
parents: 14383
diff changeset
401
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
402 %!function f = __rosenb__ (x)
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
403 %! n = length (x);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
404 %! f = sumsq (1 - x(1:n-1)) + 100 * sumsq (x(2:n) - x(1:n-1).^2);
13305
63463570d9fe Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents: 13027
diff changeset
405 %!endfunction
19833
9fc020886ae9 maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19697
diff changeset
406 %!
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
407 %!test
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
408 %! [x, fval, info, out] = fminunc (@__rosenb__, [5, -5]);
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
409 %! tol = 2e-5;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
410 %! assert (info > 0);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
411 %! assert (x, ones (1, 2), tol);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
412 %! assert (fval, 0, tol);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
413 %!test
20695
6faaab833605 fminunc.m: Clean up function to meet Octave coding standards.
Rik <rik@octave.org>
parents: 20694
diff changeset
414 %! [x, fval, info, out] = fminunc (@__rosenb__, zeros (1, 4));
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
415 %! tol = 2e-5;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
416 %! assert (info > 0);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
417 %! assert (x, ones (1, 4), tol);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
418 %! assert (fval, 0, tol);
19833
9fc020886ae9 maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19697
diff changeset
419
9fc020886ae9 maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19697
diff changeset
420 ## Test FunValCheck works correctly
14386
59aab666f2bf Extend "FunValCheck" option to optimization routines to detect Inf values.
Rik <octave@nomad.inbox5.com>
parents: 14383
diff changeset
421 %!assert (fminunc (@(x) x^2, 1, optimset ("FunValCheck", "on")), 0, eps)
59aab666f2bf Extend "FunValCheck" option to optimization routines to detect Inf values.
Rik <octave@nomad.inbox5.com>
parents: 14383
diff changeset
422 %!error <non-real value> fminunc (@(x) x + i, 1, optimset ("FunValCheck", "on"))
59aab666f2bf Extend "FunValCheck" option to optimization routines to detect Inf values.
Rik <octave@nomad.inbox5.com>
parents: 14383
diff changeset
423 %!error <NaN value> fminunc (@(x) x + NaN, 1, optimset ("FunValCheck", "on"))
59aab666f2bf Extend "FunValCheck" option to optimization routines to detect Inf values.
Rik <octave@nomad.inbox5.com>
parents: 14383
diff changeset
424 %!error <Inf value> fminunc (@(x) x + Inf, 1, optimset ("FunValCheck", "on"))
59aab666f2bf Extend "FunValCheck" option to optimization routines to detect Inf values.
Rik <octave@nomad.inbox5.com>
parents: 14383
diff changeset
425
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
426
9899
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
427 ## Solve the double dogleg trust-region minimization problem:
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
428 ## Minimize 1/2*norm(r*x)^2 subject to the constraint norm(d.*x) <= delta,
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
429 ## 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: 9758
diff changeset
430
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
431 ## TODO: error checks
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
432 ## TODO: handle singularity, or leave it up to mldivide?
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
433
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
434 function x = __doglegm__ (r, g, d, delta)
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
435 ## Get Gauss-Newton direction.
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
436 b = r' \ g;
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
437 x = r \ b;
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
438 xn = norm (d .* x);
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
439 if (xn > delta)
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
440 ## GN is too big, get scaled gradient.
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
441 s = g ./ d;
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
442 sn = norm (s);
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
443 if (sn > 0)
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
444 ## Normalize and rescale.
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
445 s = (s / sn) ./ d;
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
446 ## Get the line minimizer in s direction.
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
447 tn = norm (r*s);
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
448 snm = (sn / tn) / tn;
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
449 if (snm < delta)
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10201
diff changeset
450 ## Get the dogleg path minimizer.
9899
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
451 bn = norm (b);
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
452 dxn = delta/xn; snmd = snm/delta;
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
453 t = (bn/sn) * (bn/xn) * snmd;
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
454 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: 9758
diff changeset
455 alpha = dxn*(1-snmd^2) / t;
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
456 else
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
457 alpha = 0;
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
458 endif
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
459 else
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
460 alpha = delta / xn;
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
461 snm = 0;
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
462 endif
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
463 ## Form the appropriate convex combination.
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
464 x = alpha * x + ((1-alpha) * min (snm, delta)) * s;
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
465 endif
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
466 endfunction
17338
1c89599167a6 maint: End m-files with 1 blank line.
Rik <rik@octave.org>
parents: 17281
diff changeset
467