annotate scripts/optimization/fminunc.m @ 9212:6feb27c38da1

support central differences in fminunc and fsolve
author Jaroslav Hajek <highegg@gmail.com>
date Mon, 18 May 2009 09:46:49 +0200
parents 25f50d2d76b3
children e598248a060d
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
1 ## Copyright (C) 2008, 2009 VZLU Prague, a.s.
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 -*-
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
22 ## @deftypefn{Function File} {} fminunc (@var{fcn}, @var{x0}, @var{options})
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
23 ## @deftypefnx{Function File} {[@var{x}, @var{fvec}, @var{info}, @var{output}, @var{fjac}]} = fminunc (@var{fcn}, @dots{})
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
24 ## Solve a unconstrained optimization problem defined by the function @var{fcn}.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
25 ## @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
26 ## and return the objective function value, optionally with gradient.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
27 ## In other words, this function attempts to determine a vector @var{x} such
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
28 ## that @code{@var{fcn} (@var{x})} is a local minimum.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
29 ## @var{x0} determines a starting guess. The shape of @var{x0} is preserved
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
30 ## in all calls to @var{fcn}, but otherwise it is treated as a column vector.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
31 ## @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
32 ## Currently, @code{fminunc} recognizes these options:
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
33 ## @code{"FunValCheck"}, @code{"OutputFcn"}, @code{"TolX"},
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
34 ## @code{"TolFun"}, @code{"MaxIter"}, @code{"MaxFunEvals"},
9212
6feb27c38da1 support central differences in fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9207
diff changeset
35 ## @code{"GradObj"}, @code{"FinDiffType"}.
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
36 ##
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
37 ## 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
38 ## 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
39 ## of right-hand sides at the requested point. @code{"TolX"} specifies
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
40 ## the termination tolerance in the unknown variables, while
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
41 ## @code{"TolFun"} is a tolerance for equations. Default is @code{1e-7}
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
42 ## for both @code{"TolX"} and @code{"TolFun"}.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
43 ##
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
44 ## 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
45 ##
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
46 ## 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
47 ## evaluated at @var{x}, and @var{info} may be one of the following values:
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
48 ##
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
49 ## @table @asis
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
50 ## @item 1
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
51 ## Converged to a solution point. Relative gradient error is less than specified
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
52 ## by TolFun.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
53 ## @item 2
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
54 ## Last relative step size was less that TolX.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
55 ## @item 3
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
56 ## Last relative decrease in func value was less than TolF.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
57 ## @item 0
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
58 ## Iteration limit exceeded.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
59 ## @item -3
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
60 ## The trust region radius became excessively small.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
61 ## @end table
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
62 ##
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
63 ## Note: If you only have a single nonlinear equation of one variable, using
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
64 ## @code{fminbnd} is usually a much better idea.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
65 ## @seealso{fminbnd, optimset}
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
66 ## @end deftypefn
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
67
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
68 ## PKG_ADD: __all_opts__ ("fminunc");
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
69
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
70 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
71
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
72 ## Get default options if requested.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
73 if (nargin == 1 && ischar (fcn) && strcmp (fcn, 'defaults'))
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
74 x = optimset ("MaxIter", 400, "MaxFunEvals", Inf, \
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
75 "GradObj", "off", "TolX", 1.5e-8, "TolFun", 1.5e-8,
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
76 "OutputFcn", [], "FunValCheck", "off",
9212
6feb27c38da1 support central differences in fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9207
diff changeset
77 "FinDiffType", "central");
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
78 return;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
79 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
80
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
81 if (nargin < 2 || nargin > 3 || ! ismatrix (x0))
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
82 print_usage ();
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
83 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
84
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
85 if (ischar (fcn))
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
86 fcn = str2func (fcn);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
87 endif
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 xsiz = size (x0);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
90 n = numel (x0);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
91
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
92 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
93 cdif = strcmpi (optimget (options, "FinDiffType", "central"), "central");
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
94 maxiter = optimget (options, "MaxIter", 400);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
95 maxfev = optimget (options, "MaxFunEvals", Inf);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
96 outfcn = optimget (options, "OutputFcn");
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 funvalchk = strcmpi (optimget (options, "FunValCheck", "off"), "on");
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
99
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
100 if (funvalchk)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
101 ## Replace fcn with a guarded version.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
102 fcn = @(x) guarded_eval (fcn, x);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
103 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
104
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
105 ## 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
106 ## prefers accuracy to performance.
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 macheps = eps (class (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 tolx = optimget (options, "TolX", sqrt (macheps));
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
111 tolf = optimget (options, "TolFun", sqrt (macheps));
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
112
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
113 factor = 100;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
114 ## FIXME: TypicalX corresponds to user scaling (???)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
115 autodg = true;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
116
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
117 niter = 1;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
118 nfev = 0;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
119
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
120 x = x0(:);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
121 info = 0;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
122
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
123 ## Initial evaluation.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
124 fval = fcn (reshape (x, xsiz));
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
125 n = length (x);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
126
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
127 if (! isempty (outfcn))
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
128 optimvalues.iter = niter;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
129 optimvalues.funccount = nfev;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
130 optimvalues.fval = fval;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
131 optimvalues.searchdirection = zeros (n, 1);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
132 state = 'init';
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
133 stop = outfcn (x, optimvalues, state);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
134 if (stop)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
135 info = -1;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
136 break;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
137 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
138 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
139
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
140 nsuciter = 0;
9199
399884c9d4a1 import the step adaptation strategy from fsolve to fminunc
Jaroslav Hajek <highegg@gmail.com>
parents: 9084
diff changeset
141 lastratio = 0;
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
142
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
143 grad = [];
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
144
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
145 ## Outer loop.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
146 while (niter < maxiter && nfev < maxfev && ! info)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
147
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
148 grad0 = grad;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
149
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
150 ## Calculate function value and gradient (possibly via FD).
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
151 if (has_grad)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
152 [fval, grad] = fcn (reshape (x, xsiz));
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
153 grad = grad(:);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
154 nfev ++;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
155 else
9212
6feb27c38da1 support central differences in fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9207
diff changeset
156 grad = __fdjac__ (fcn, reshape (x, xsiz), fval, cdif)(:);
6feb27c38da1 support central differences in fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9207
diff changeset
157 nfev += (1 + cdif) * length (x);
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
158 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
159
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
160 if (niter == 1)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
161 ## Initialize by identity matrix.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
162 hesr = eye (n);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
163 else
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
164 ## Use the damped BFGS formula.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
165 y = grad - grad0;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
166 sBs = sumsq (w);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
167 Bs = hesr'*w;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
168 sy = y'*s;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
169 if (sy >= 0.2*sBs)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
170 theta = 1;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
171 else
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
172 theta = 0.8*sBs / (sBs - sy);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
173 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
174 r = theta * y + (1-theta) * Bs;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
175 hesr = cholupdate (hesr, r / sqrt (s'*r), "+");
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
176 [hesr, info] = cholupdate (hesr, Bs / sqrt (sBs), "-");
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
177 if (info)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
178 hesr = eye (n);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
179 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
180 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
181
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
182 ## Second derivatives approximate the hessian.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
183 d2f = norm (hesr, 'columns').';
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
184 if (niter == 1)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
185 dg = d2f;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
186 xn = norm (dg .* x);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
187 ## FIXME: something better?
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
188 delta = max (factor * xn, 1);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
189 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
190
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
191 ## FIXME: maybe fixed lower and upper bounds?
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
192 dg = max (0.1*dg, d2f);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
193
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
194 ## 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
195 ## 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
196 ## 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
197 if (norm (grad) <= tolf*n*xn)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
198 info = 1;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
199 break;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
200 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
201
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
202 suc = false;
9207
25f50d2d76b3 improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9199
diff changeset
203 decfac = 0.5;
25f50d2d76b3 improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9199
diff changeset
204
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
205 ## Inner loop.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
206 while (! suc && niter <= maxiter && nfev < maxfev && ! info)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
207
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
208 s = - __dogleg__ (hesr, grad, dg, delta, true);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
209
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
210 sn = norm (dg .* s);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
211 if (niter == 1)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
212 delta = min (delta, sn);
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
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
215 fval1 = fcn (reshape (x + s, xsiz)) (:);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
216 nfev ++;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
217
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
218 if (fval1 < fval)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
219 ## Scaled actual reduction.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
220 actred = (fval - fval1) / (abs (fval1) + abs (fval));
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
221 else
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
222 actred = -1;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
223 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
224
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
225 w = hesr*s;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
226 ## Scaled predicted reduction, and ratio.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
227 t = 1/2 * sumsq (w) + grad'*s;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
228 if (t < 0)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
229 prered = -t/(abs (fval) + abs (fval + t));
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
230 ratio = actred / prered;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
231 else
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
232 prered = 0;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
233 ratio = 0;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
234 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
235
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
236 ## Update delta.
9207
25f50d2d76b3 improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9199
diff changeset
237 if (ratio < min(max(0.1, 0.8*lastratio), 0.9))
25f50d2d76b3 improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9199
diff changeset
238 delta *= decfac;
25f50d2d76b3 improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9199
diff changeset
239 decfac ^= 1.4142;
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
240 if (delta <= 1e1*macheps*xn)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
241 ## Trust region became uselessly small.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
242 info = -3;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
243 break;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
244 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
245 else
9199
399884c9d4a1 import the step adaptation strategy from fsolve to fminunc
Jaroslav Hajek <highegg@gmail.com>
parents: 9084
diff changeset
246 lastratio = ratio;
9207
25f50d2d76b3 improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9199
diff changeset
247 decfac = 0.5;
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
248 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
249 delta = 1.4142*sn;
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
250 elseif (ratio >= 0.5)
9207
25f50d2d76b3 improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9199
diff changeset
251 delta = max (delta, 1.4142*sn);
9084
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 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 if (ratio >= 1e-4)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
256 ## Successful iteration.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
257 x += s;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
258 xn = norm (dg .* x);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
259 fval = fval1;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
260 nsuciter ++;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
261 suc = true;
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
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
264 niter ++;
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 ## 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
267 if (! isempty (outfcn))
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
268 optimvalues.iter = niter;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
269 optimvalues.funccount = nfev;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
270 optimvalues.fval = fval;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
271 optimvalues.searchdirection = s;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
272 state = 'iter';
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
273 stop = outfcn (x, optimvalues, state);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
274 if (stop)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
275 info = -1;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
276 break;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
277 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
278 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
279
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
280 ## Tests for termination conditions. A mysterious place, anything
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
281 ## can happen if you change something here...
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
282
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
283 ## 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
284 ## 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
285 ## 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
286 ## iterations, so we need scaling-independent tolerances wherever
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
287 ## possible.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
288
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
289 ## The following tests done only after successful step.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
290 if (ratio >= 1e-4)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
291 ## 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
292 ## but compare to scaled step, so nothing bad.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
293 if (sn <= tolx*xn)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
294 info = 2;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
295 ## Again a classic one.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
296 elseif (actred < tolf)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
297 info = 3;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
298 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
299 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
300
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
301 endwhile
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
302 endwhile
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
303
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
304 ## Restore original shapes.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
305 x = reshape (x, xsiz);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
306
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
307 output.iterations = niter;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
308 output.successful = nsuciter;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
309 output.funcCount = nfev;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
310
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
311 endfunction
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
312
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
313 ## 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
314 ## bad results.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
315 function [fx, gx] = guarded_eval (fun, x)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
316 if (nargout > 1)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
317 [fx, gx] = fun (x);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
318 else
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
319 fx = fun (x);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
320 gx = [];
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
321 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
322
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
323 if (! (isreal (fx) && isreal (jx)))
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
324 error ("fminunc:notreal", "fminunc: non-real value encountered");
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
325 elseif (complexeqn && ! (isnumeric (fx) && isnumeric(jx)))
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
326 error ("fminunc:notnum", "fminunc: non-numeric value encountered");
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
327 elseif (any (isnan (fx(:))))
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
328 error ("fminunc:isnan", "fminunc: NaN value encountered");
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 endfunction
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
331
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
332 %!function f = rosenb (x)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
333 %! n = length (x);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
334 %! f = sumsq (1 - x(1:n-1)) + 100 * sumsq (x(2:n) - x(1:n-1).^2);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
335 %!test
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
336 %! [x, fval, info, out] = fminunc (@rosenb, [5, -5]);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
337 %! tol = 2e-5;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
338 %! assert (info > 0);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
339 %! assert (x, ones (1, 2), tol);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
340 %! assert (fval, 0, tol);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
341 %!test
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
342 %! [x, fval, info, out] = fminunc (@rosenb, zeros (1, 4));
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
343 %! tol = 2e-5;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
344 %! assert (info > 0);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
345 %! assert (x, ones (1, 4), tol);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
346 %! assert (fval, 0, tol);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
347