Mercurial > hg > octave-nkf
annotate scripts/optimization/fminunc.m @ 15107:03381a36f70d
generate convenience libraries for new parse-tree and interpfcn subdirectories
* src/Makefile.am (liboctinterp_la_SOURCES): Include octave.cc in the
list, not $(DIST_SRC).
(liboctinterp_la_LIBADD): Include octave-value/liboctave-value.la,
parse-tree/libparse-tree.la, interp-core/libinterp-core.la,
interpfcn/libinterpfcn.la, and corefcn/libcorefcn.la in the list.
* src/interp-core/module.mk (noinst_LTLIBRARIES): Add
interp-core/libinterp-core.la to the list.
(interp_core_libinterp_core_la_SOURCES): New variable.
* src/interpfcn/module.mk (noinst_LTLIBRARIES): Add
interpfcn/libinterpfcn.la to the list.
(interpfcn_libinterpfcn_la_SOURCES): New variable.
* src/parse-tree/module.mk (noinst_LTLIBRARIES): Add
parse-tree/libparse-tree.la to the list.
(parse_tree_libparse_tree_la_SOURCES): New variable.
* src/octave-value/module.mk (noinst_LTLIBRARIES): Add
octave-value/liboctave-value.la to the list.
(octave_value_liboctave_value_la_SOURCES): New variable.
author | John W. Eaton <jwe@octave.org> |
---|---|
date | Sun, 05 Aug 2012 09:04:30 -0400 |
parents | e0525ecf156e |
children | bc924baa2c4e |
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. |
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}) |
12578
f5a780d675a1
Clean up operator and function indices in documentation.
Rik <octave@nomad.inbox5.com>
parents:
12576
diff
changeset
|
24 ## @deftypefnx {Function File} {[@var{x}, @var{fvec}, @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 ## |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
28 ## @var{fcn} should accepts a vector (array) defining the unknown variables, |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
29 ## and return the objective function value, optionally with gradient. |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
30 ## In other words, this function attempts to determine a vector @var{x} such |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
31 ## that @code{@var{fcn} (@var{x})} is a local minimum. |
9758
09da0bd91412
Periodic grammar check of Octave documentation files to ensure common format
Rik <rdrider0-list@yahoo.com>
parents:
9633
diff
changeset
|
32 ## @var{x0} determines a starting guess. The shape of @var{x0} is preserved |
10687
a8ce6bdecce5
Improve documentation strings.
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
33 ## in all calls to @var{fcn}, but otherwise is treated as a column vector. |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
34 ## @var{options} is a structure specifying additional options. |
9212
6feb27c38da1
support central differences in fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9207
diff
changeset
|
35 ## Currently, @code{fminunc} recognizes these options: |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
36 ## @code{"FunValCheck"}, @code{"OutputFcn"}, @code{"TolX"}, |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
37 ## @code{"TolFun"}, @code{"MaxIter"}, @code{"MaxFunEvals"}, |
10201
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
38 ## @code{"GradObj"}, @code{"FinDiffType"}, |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
39 ## @code{"TypicalX"}, @code{"AutoScaling"}. |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
40 ## |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
41 ## If @code{"GradObj"} is @code{"on"}, it specifies that @var{fcn}, |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
42 ## called with 2 output arguments, also returns the Jacobian matrix |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
43 ## of right-hand sides at the requested point. @code{"TolX"} specifies |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
44 ## the termination tolerance in the unknown variables, while |
9758
09da0bd91412
Periodic grammar check of Octave documentation files to ensure common format
Rik <rdrider0-list@yahoo.com>
parents:
9633
diff
changeset
|
45 ## @code{"TolFun"} is a tolerance for equations. Default is @code{1e-7} |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
46 ## for both @code{"TolX"} and @code{"TolFun"}. |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
47 ## |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
48 ## For description of the other options, see @code{optimset}. |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
49 ## |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
50 ## On return, @var{fval} contains the value of the function @var{fcn} |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
51 ## 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
|
52 ## |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
53 ## @table @asis |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
54 ## @item 1 |
10793
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
10687
diff
changeset
|
55 ## Converged to a solution point. Relative gradient error is less than |
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
10687
diff
changeset
|
56 ## specified |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
57 ## by TolFun. |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
58 ## |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
59 ## @item 2 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
60 ## Last relative step size was less that TolX. |
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 3 |
12576
a1e386b9ef4b
Spellcheck documentation for 3.4.1 release.
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
63 ## Last relative decrease in function value was less than TolF. |
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 0 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
66 ## Iteration limit exceeded. |
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 -3 |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
69 ## The trust region radius became excessively small. |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
70 ## @end table |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
71 ## |
9627
5bcfa0b346e8
fix extra outputs in fminunc
Jaroslav Hajek <highegg@gmail.com>
parents:
9623
diff
changeset
|
72 ## Optionally, fminunc can also yield a structure with convergence statistics |
10793
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
10687
diff
changeset
|
73 ## (@var{output}), the output gradient (@var{grad}) and approximate Hessian |
9627
5bcfa0b346e8
fix extra outputs in fminunc
Jaroslav Hajek <highegg@gmail.com>
parents:
9623
diff
changeset
|
74 ## (@var{hess}). |
5bcfa0b346e8
fix extra outputs in fminunc
Jaroslav Hajek <highegg@gmail.com>
parents:
9623
diff
changeset
|
75 ## |
14895
e0525ecf156e
Add new function fminsearch.m
Andy Adler <andy@analyti.ca>
parents:
14868
diff
changeset
|
76 ## Notes: If you only have a single nonlinear equation of one variable then |
e0525ecf156e
Add new function fminsearch.m
Andy Adler <andy@analyti.ca>
parents:
14868
diff
changeset
|
77 ## using @code{fminbnd} is usually a much better idea. The algorithm used is a |
e0525ecf156e
Add new function fminsearch.m
Andy Adler <andy@analyti.ca>
parents:
14868
diff
changeset
|
78 ## gradient search which depends on the objective function being differentiable. |
e0525ecf156e
Add new function fminsearch.m
Andy Adler <andy@analyti.ca>
parents:
14868
diff
changeset
|
79 ## If the function has discontinuities it may be better to use a derivative-free |
e0525ecf156e
Add new function fminsearch.m
Andy Adler <andy@analyti.ca>
parents:
14868
diff
changeset
|
80 ## algorithm such as @code{fminsearch}. |
e0525ecf156e
Add new function fminsearch.m
Andy Adler <andy@analyti.ca>
parents:
14868
diff
changeset
|
81 ## @seealso{fminbnd, fminsearch, optimset} |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
82 ## @end deftypefn |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
83 |
13027
b9a89ca0fb75
prevent optimization functions from setting ans in workspace at startup
John W. Eaton <jwe@octave.org>
parents:
12578
diff
changeset
|
84 ## 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
|
85 ## PKG_ADD: [~] = __all_opts__ ("fminunc"); |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
86 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
87 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
|
88 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
89 ## Get default options if requested. |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
90 if (nargin == 1 && ischar (fcn) && strcmp (fcn, 'defaults')) |
14552
86854d032a37
maint: miscellaneous style fixes for .m files
John W. Eaton <jwe@octave.org>
parents:
14386
diff
changeset
|
91 x = optimset ("MaxIter", 400, "MaxFunEvals", Inf, |
86854d032a37
maint: miscellaneous style fixes for .m files
John W. Eaton <jwe@octave.org>
parents:
14386
diff
changeset
|
92 "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
|
93 "OutputFcn", [], "FunValCheck", "off", |
86854d032a37
maint: miscellaneous style fixes for .m files
John W. Eaton <jwe@octave.org>
parents:
14386
diff
changeset
|
94 "FinDiffType", "central", |
86854d032a37
maint: miscellaneous style fixes for .m files
John W. Eaton <jwe@octave.org>
parents:
14386
diff
changeset
|
95 "TypicalX", [], "AutoScaling", "off"); |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
96 return; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
97 endif |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
98 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
99 if (nargin < 2 || nargin > 3 || ! ismatrix (x0)) |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
100 print_usage (); |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
101 endif |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
102 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
103 if (ischar (fcn)) |
9464
e598248a060d
safer str2func use in optim functions
Jaroslav Hajek <highegg@gmail.com>
parents:
9212
diff
changeset
|
104 fcn = str2func (fcn, "global"); |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
105 endif |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
106 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
107 xsiz = size (x0); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
108 n = numel (x0); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
109 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
110 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
|
111 cdif = strcmpi (optimget (options, "FinDiffType", "central"), "central"); |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
112 maxiter = optimget (options, "MaxIter", 400); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
113 maxfev = optimget (options, "MaxFunEvals", Inf); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
114 outfcn = optimget (options, "OutputFcn"); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
115 |
10201
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
116 ## Get scaling matrix using the TypicalX option. If set to "auto", the |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
117 ## scaling matrix is estimated using the jacobian. |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
118 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
|
119 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
|
120 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
|
121 endif |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
122 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
|
123 if (! autoscale) |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
124 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
|
125 endif |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
126 |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
127 funvalchk = strcmpi (optimget (options, "FunValCheck", "off"), "on"); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
128 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
129 if (funvalchk) |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
130 ## Replace fcn with a guarded version. |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
131 fcn = @(x) guarded_eval (fcn, x); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
132 endif |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
133 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
134 ## These defaults are rather stringent. I think that normally, user |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
135 ## prefers accuracy to performance. |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
136 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
137 macheps = eps (class (x0)); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
138 |
10201
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
139 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
|
140 tolf = optimget (options, "TolFun", 1e-7); |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
141 |
9623
bc0739d02724
update initial TR step for fsolve and fminunc
Jaroslav Hajek <highegg@gmail.com>
parents:
9464
diff
changeset
|
142 factor = 0.1; |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
143 ## FIXME: TypicalX corresponds to user scaling (???) |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
144 autodg = true; |
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 niter = 1; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
147 nfev = 0; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
148 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
149 x = x0(:); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
150 info = 0; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
151 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
152 ## Initial evaluation. |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
153 fval = fcn (reshape (x, xsiz)); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
154 n = length (x); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
155 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
156 if (! isempty (outfcn)) |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
157 optimvalues.iter = niter; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
158 optimvalues.funccount = nfev; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
159 optimvalues.fval = fval; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
160 optimvalues.searchdirection = zeros (n, 1); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
161 state = 'init'; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
162 stop = outfcn (x, optimvalues, state); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
163 if (stop) |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
164 info = -1; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
165 break; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
166 endif |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
167 endif |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
168 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
169 nsuciter = 0; |
9199
399884c9d4a1
import the step adaptation strategy from fsolve to fminunc
Jaroslav Hajek <highegg@gmail.com>
parents:
9084
diff
changeset
|
170 lastratio = 0; |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
171 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
172 grad = []; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
173 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
174 ## Outer loop. |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
175 while (niter < maxiter && nfev < maxfev && ! info) |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
176 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
177 grad0 = grad; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
178 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
179 ## Calculate function value and gradient (possibly via FD). |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
180 if (has_grad) |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
181 [fval, grad] = fcn (reshape (x, xsiz)); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
182 grad = grad(:); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
183 nfev ++; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
184 else |
10201
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
185 grad = __fdjac__ (fcn, reshape (x, xsiz), fval, typicalx, cdif)(:); |
9212
6feb27c38da1
support central differences in fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9207
diff
changeset
|
186 nfev += (1 + cdif) * length (x); |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
187 endif |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
188 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
189 if (niter == 1) |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
190 ## Initialize by identity matrix. |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
191 hesr = eye (n); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
192 else |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
193 ## Use the damped BFGS formula. |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
194 y = grad - grad0; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
195 sBs = sumsq (w); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
196 Bs = hesr'*w; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
197 sy = y'*s; |
9633
ecc2c556f844
simplify damped BFGS formula
Jaroslav Hajek <highegg@gmail.com>
parents:
9631
diff
changeset
|
198 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
|
199 r = theta * y + (1-theta) * Bs; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
200 hesr = cholupdate (hesr, r / sqrt (s'*r), "+"); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
201 [hesr, info] = cholupdate (hesr, Bs / sqrt (sBs), "-"); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
202 if (info) |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
203 hesr = eye (n); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
204 endif |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
205 endif |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
206 |
10201
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
207 if (autoscale) |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
208 ## Second derivatives approximate the hessian. |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
209 d2f = norm (hesr, 'columns').'; |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
210 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
|
211 dg = d2f; |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
212 else |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
213 ## 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
|
214 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
|
215 endif |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
216 endif |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
217 |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
218 if (niter == 1) |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
219 xn = norm (dg .* x); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
220 ## FIXME: something better? |
9628
73e6ad869f08
further correct initial TR step strategy
Jaroslav Hajek <highegg@gmail.com>
parents:
9627
diff
changeset
|
221 delta = factor * max (xn, 1); |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
222 endif |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
223 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
224 ## FIXME -- why tolf*n*xn? If abs (e) ~ abs(x) * eps is a vector |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
225 ## 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
|
226 ## 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
|
227 if (norm (grad) <= tolf*n*xn) |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
228 info = 1; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
229 break; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
230 endif |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
231 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
232 suc = false; |
9207
25f50d2d76b3
improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9199
diff
changeset
|
233 decfac = 0.5; |
25f50d2d76b3
improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9199
diff
changeset
|
234 |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
235 ## Inner loop. |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
236 while (! suc && niter <= maxiter && nfev < maxfev && ! info) |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
237 |
9631
00958d0c4e3c
split __dogleg__ > __doglegm__
Jaroslav Hajek <highegg@gmail.com>
parents:
9628
diff
changeset
|
238 s = - __doglegm__ (hesr, grad, dg, delta); |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
239 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
240 sn = norm (dg .* s); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
241 if (niter == 1) |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
242 delta = min (delta, sn); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
243 endif |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
244 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
245 fval1 = fcn (reshape (x + s, xsiz)) (:); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
246 nfev ++; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
247 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
248 if (fval1 < fval) |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
249 ## Scaled actual reduction. |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
250 actred = (fval - fval1) / (abs (fval1) + abs (fval)); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
251 else |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
252 actred = -1; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
253 endif |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
254 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
255 w = hesr*s; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
256 ## Scaled predicted reduction, and ratio. |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
257 t = 1/2 * sumsq (w) + grad'*s; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
258 if (t < 0) |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
259 prered = -t/(abs (fval) + abs (fval + t)); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
260 ratio = actred / prered; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
261 else |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
262 prered = 0; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
263 ratio = 0; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
264 endif |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
265 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
266 ## Update delta. |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14552
diff
changeset
|
267 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
|
268 delta *= decfac; |
25f50d2d76b3
improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9199
diff
changeset
|
269 decfac ^= 1.4142; |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
270 if (delta <= 1e1*macheps*xn) |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
271 ## Trust region became uselessly small. |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
272 info = -3; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
273 break; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
274 endif |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
275 else |
9199
399884c9d4a1
import the step adaptation strategy from fsolve to fminunc
Jaroslav Hajek <highegg@gmail.com>
parents:
9084
diff
changeset
|
276 lastratio = ratio; |
9207
25f50d2d76b3
improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9199
diff
changeset
|
277 decfac = 0.5; |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
278 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
|
279 delta = 1.4142*sn; |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
280 elseif (ratio >= 0.5) |
9207
25f50d2d76b3
improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9199
diff
changeset
|
281 delta = max (delta, 1.4142*sn); |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
282 endif |
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 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
285 if (ratio >= 1e-4) |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
286 ## Successful iteration. |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
287 x += s; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
288 xn = norm (dg .* x); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
289 fval = fval1; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
290 nsuciter ++; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
291 suc = true; |
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 niter ++; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
295 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
296 ## FIXME: should outputfcn be only called after a successful iteration? |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
297 if (! isempty (outfcn)) |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
298 optimvalues.iter = niter; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
299 optimvalues.funccount = nfev; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
300 optimvalues.fval = fval; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
301 optimvalues.searchdirection = s; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
302 state = 'iter'; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
303 stop = outfcn (x, optimvalues, state); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
304 if (stop) |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
305 info = -1; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
306 break; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
307 endif |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
308 endif |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
309 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
310 ## Tests for termination conditions. A mysterious place, anything |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
311 ## 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
|
312 |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
313 ## 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
|
314 ## is that for a tolerance that depends on scaling, only 0 makes |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
315 ## sense as a default value. But 0 usually means uselessly long |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
316 ## iterations, so we need scaling-independent tolerances wherever |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
317 ## possible. |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
318 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
319 ## The following tests done only after successful step. |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
320 if (ratio >= 1e-4) |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
321 ## This one is classic. Note that we use scaled variables again, |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
322 ## but compare to scaled step, so nothing bad. |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
323 if (sn <= tolx*xn) |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
324 info = 2; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
325 ## Again a classic one. |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
326 elseif (actred < tolf) |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
327 info = 3; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
328 endif |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
329 endif |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
330 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
331 endwhile |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
332 endwhile |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
333 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
334 ## Restore original shapes. |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
335 x = reshape (x, xsiz); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
336 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
337 output.iterations = niter; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
338 output.successful = nsuciter; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
339 output.funcCount = nfev; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
340 |
9627
5bcfa0b346e8
fix extra outputs in fminunc
Jaroslav Hajek <highegg@gmail.com>
parents:
9623
diff
changeset
|
341 if (nargout > 5) |
5bcfa0b346e8
fix extra outputs in fminunc
Jaroslav Hajek <highegg@gmail.com>
parents:
9623
diff
changeset
|
342 hess = hesr'*hesr; |
5bcfa0b346e8
fix extra outputs in fminunc
Jaroslav Hajek <highegg@gmail.com>
parents:
9623
diff
changeset
|
343 endif |
5bcfa0b346e8
fix extra outputs in fminunc
Jaroslav Hajek <highegg@gmail.com>
parents:
9623
diff
changeset
|
344 |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
345 endfunction |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
346 |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
347 ## An assistant function that evaluates a function handle and checks for |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
348 ## bad results. |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
349 function [fx, gx] = guarded_eval (fun, x) |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
350 if (nargout > 1) |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
351 [fx, gx] = fun (x); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
352 else |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
353 fx = fun (x); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
354 gx = []; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
355 endif |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
356 |
14383
07c55bceca23
Fix guarded_eval() subfunction in fminunc (bug #35534).
Olaf Till <olaf.till@uni-jena.de>
parents:
14138
diff
changeset
|
357 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
|
358 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
|
359 elseif (any (isnan (fx(:)))) |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
360 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
|
361 elseif (any (isinf (fx(:)))) |
59aab666f2bf
Extend "FunValCheck" option to optimization routines to detect Inf values.
Rik <octave@nomad.inbox5.com>
parents:
14383
diff
changeset
|
362 error ("fminunc:isinf", "fminunc: Inf value encountered"); |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
363 endif |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
364 endfunction |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
365 |
14386
59aab666f2bf
Extend "FunValCheck" option to optimization routines to detect Inf values.
Rik <octave@nomad.inbox5.com>
parents:
14383
diff
changeset
|
366 |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
367 %!function f = __rosenb (x) |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
368 %! n = length (x); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
369 %! 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
|
370 %!endfunction |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
371 %!test |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
372 %! [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
|
373 %! tol = 2e-5; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
374 %! assert (info > 0); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
375 %! assert (x, ones (1, 2), tol); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
376 %! assert (fval, 0, tol); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
377 %!test |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
378 %! [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
|
379 %! tol = 2e-5; |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
380 %! assert (info > 0); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
381 %! assert (x, ones (1, 4), tol); |
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
382 %! assert (fval, 0, tol); |
14386
59aab666f2bf
Extend "FunValCheck" option to optimization routines to detect Inf values.
Rik <octave@nomad.inbox5.com>
parents:
14383
diff
changeset
|
383 %% Test FunValCheck works correctly |
59aab666f2bf
Extend "FunValCheck" option to optimization routines to detect Inf values.
Rik <octave@nomad.inbox5.com>
parents:
14383
diff
changeset
|
384 %!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
|
385 %!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
|
386 %!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
|
387 %!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
|
388 |
9084
b7210faa3ed0
implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
389 |
9899
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
390 ## 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
|
391 ## 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
|
392 ## 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
|
393 |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
394 ## TODO: error checks |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
395 ## 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
|
396 |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
397 function x = __doglegm__ (r, g, d, delta) |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
398 ## Get Gauss-Newton direction. |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
399 b = r' \ g; |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
400 x = r \ b; |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
401 xn = norm (d .* x); |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
402 if (xn > delta) |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
403 ## GN is too big, get scaled gradient. |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
404 s = g ./ d; |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
405 sn = norm (s); |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
406 if (sn > 0) |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
407 ## Normalize and rescale. |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
408 s = (s / sn) ./ d; |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
409 ## Get the line minimizer in s direction. |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
410 tn = norm (r*s); |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
411 snm = (sn / tn) / tn; |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
412 if (snm < delta) |
10549 | 413 ## Get the dogleg path minimizer. |
9899
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
414 bn = norm (b); |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
415 dxn = delta/xn; snmd = snm/delta; |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
416 t = (bn/sn) * (bn/xn) * snmd; |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
417 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
|
418 alpha = dxn*(1-snmd^2) / t; |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
419 else |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
420 alpha = 0; |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
421 endif |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
422 else |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
423 alpha = delta / xn; |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
424 snm = 0; |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
425 endif |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
426 ## Form the appropriate convex combination. |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
427 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
|
428 endif |
9f25290a35e8
more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
429 endfunction |