annotate scripts/optimization/fminunc.m @ 14138:72c96de7a403 stable

maint: update copyright notices for 2012
author John W. Eaton <jwe@octave.org>
date Mon, 02 Jan 2012 14:25:41 -0500
parents 63463570d9fe
children 07c55bceca23
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
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}.
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
27 ## @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
28 ## 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
29 ## 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
30 ## 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
31 ## @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
32 ## 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
33 ## @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
34 ## Currently, @code{fminunc} recognizes these options:
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
35 ## @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
36 ## @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
37 ## @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
38 ## @code{"TypicalX"}, @code{"AutoScaling"}.
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
39 ##
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
40 ## 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
41 ## 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
42 ## 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
43 ## 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
44 ## @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
45 ## 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
46 ##
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
47 ## 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
48 ##
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
49 ## 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
50 ## 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
51 ##
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
52 ## @table @asis
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
53 ## @item 1
10793
be55736a0783 Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents: 10687
diff changeset
54 ## 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
55 ## specified
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
56 ## by TolFun.
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10793
diff changeset
57 ##
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
58 ## @item 2
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
59 ## Last relative step size was less that TolX.
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10793
diff changeset
60 ##
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
61 ## @item 3
12576
a1e386b9ef4b Spellcheck documentation for 3.4.1 release.
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
62 ## 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
63 ##
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
64 ## @item 0
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
65 ## Iteration limit exceeded.
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10793
diff changeset
66 ##
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
67 ## @item -3
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
68 ## The trust region radius became excessively small.
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
69 ## @end table
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
70 ##
9627
5bcfa0b346e8 fix extra outputs in fminunc
Jaroslav Hajek <highegg@gmail.com>
parents: 9623
diff changeset
71 ## 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
72 ## (@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
73 ## (@var{hess}).
5bcfa0b346e8 fix extra outputs in fminunc
Jaroslav Hajek <highegg@gmail.com>
parents: 9623
diff changeset
74 ##
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
75 ## Note: If you only have a single nonlinear equation of one variable, using
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
76 ## @code{fminbnd} is usually a much better idea.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
77 ## @seealso{fminbnd, optimset}
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
78 ## @end deftypefn
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
79
13027
b9a89ca0fb75 prevent optimization functions from setting ans in workspace at startup
John W. Eaton <jwe@octave.org>
parents: 12578
diff changeset
80 ## 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
81 ## PKG_ADD: [~] = __all_opts__ ("fminunc");
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
82
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
83 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
84
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
85 ## Get default options if requested.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
86 if (nargin == 1 && ischar (fcn) && strcmp (fcn, 'defaults'))
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
87 x = optimset ("MaxIter", 400, "MaxFunEvals", Inf, \
10201
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
88 "GradObj", "off", "TolX", 1e-7, "TolFun", 1e-7,
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
89 "OutputFcn", [], "FunValCheck", "off",
10201
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
90 "FinDiffType", "central",
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
91 "TypicalX", [], "AutoScaling", "off");
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
92 return;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
93 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
94
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
95 if (nargin < 2 || nargin > 3 || ! ismatrix (x0))
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
96 print_usage ();
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
97 endif
9084
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 (ischar (fcn))
9464
e598248a060d safer str2func use in optim functions
Jaroslav Hajek <highegg@gmail.com>
parents: 9212
diff changeset
100 fcn = str2func (fcn, "global");
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
101 endif
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 xsiz = size (x0);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
104 n = numel (x0);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
105
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
106 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
107 cdif = strcmpi (optimget (options, "FinDiffType", "central"), "central");
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
108 maxiter = optimget (options, "MaxIter", 400);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
109 maxfev = optimget (options, "MaxFunEvals", Inf);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
110 outfcn = optimget (options, "OutputFcn");
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
111
10201
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
112 ## 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
113 ## 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
114 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
115 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
116 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
117 endif
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
118 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
119 if (! autoscale)
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
120 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
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
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
123 funvalchk = strcmpi (optimget (options, "FunValCheck", "off"), "on");
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
124
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
125 if (funvalchk)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
126 ## Replace fcn with a guarded version.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
127 fcn = @(x) guarded_eval (fcn, x);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
128 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
129
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
130 ## 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
131 ## prefers accuracy to performance.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
132
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
133 macheps = eps (class (x0));
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
134
10201
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
135 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
136 tolf = optimget (options, "TolFun", 1e-7);
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
137
9623
bc0739d02724 update initial TR step for fsolve and fminunc
Jaroslav Hajek <highegg@gmail.com>
parents: 9464
diff changeset
138 factor = 0.1;
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
139 ## FIXME: TypicalX corresponds to user scaling (???)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
140 autodg = true;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
141
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
142 niter = 1;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
143 nfev = 0;
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 x = x0(:);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
146 info = 0;
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 ## Initial evaluation.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
149 fval = fcn (reshape (x, xsiz));
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
150 n = length (x);
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 if (! isempty (outfcn))
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
153 optimvalues.iter = niter;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
154 optimvalues.funccount = nfev;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
155 optimvalues.fval = fval;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
156 optimvalues.searchdirection = zeros (n, 1);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
157 state = 'init';
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
158 stop = outfcn (x, optimvalues, state);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
159 if (stop)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
160 info = -1;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
161 break;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
162 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
163 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
164
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
165 nsuciter = 0;
9199
399884c9d4a1 import the step adaptation strategy from fsolve to fminunc
Jaroslav Hajek <highegg@gmail.com>
parents: 9084
diff changeset
166 lastratio = 0;
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
167
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
168 grad = [];
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
169
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
170 ## Outer loop.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
171 while (niter < maxiter && nfev < maxfev && ! info)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
172
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
173 grad0 = grad;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
174
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
175 ## Calculate function value and gradient (possibly via FD).
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
176 if (has_grad)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
177 [fval, grad] = fcn (reshape (x, xsiz));
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
178 grad = grad(:);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
179 nfev ++;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
180 else
10201
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
181 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
182 nfev += (1 + cdif) * length (x);
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
183 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
184
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
185 if (niter == 1)
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
186 ## Initialize by identity matrix.
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
187 hesr = eye (n);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
188 else
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
189 ## Use the damped BFGS formula.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
190 y = grad - grad0;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
191 sBs = sumsq (w);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
192 Bs = hesr'*w;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
193 sy = y'*s;
9633
ecc2c556f844 simplify damped BFGS formula
Jaroslav Hajek <highegg@gmail.com>
parents: 9631
diff changeset
194 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
195 r = theta * y + (1-theta) * Bs;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
196 hesr = cholupdate (hesr, r / sqrt (s'*r), "+");
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
197 [hesr, info] = cholupdate (hesr, Bs / sqrt (sBs), "-");
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
198 if (info)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
199 hesr = eye (n);
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 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
202
10201
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
203 if (autoscale)
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
204 ## 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
205 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
206 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
207 dg = d2f;
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
208 else
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
209 ## 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
210 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
211 endif
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
212 endif
5c66978f3fdf support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents: 9899
diff changeset
213
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
214 if (niter == 1)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
215 xn = norm (dg .* x);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
216 ## FIXME: something better?
9628
73e6ad869f08 further correct initial TR step strategy
Jaroslav Hajek <highegg@gmail.com>
parents: 9627
diff changeset
217 delta = factor * max (xn, 1);
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
218 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
219
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
220 ## 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
221 ## 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
222 ## 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
223 if (norm (grad) <= tolf*n*xn)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
224 info = 1;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
225 break;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
226 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
227
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
228 suc = false;
9207
25f50d2d76b3 improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9199
diff changeset
229 decfac = 0.5;
25f50d2d76b3 improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9199
diff changeset
230
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
231 ## Inner loop.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
232 while (! suc && niter <= maxiter && nfev < maxfev && ! info)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
233
9631
00958d0c4e3c split __dogleg__ > __doglegm__
Jaroslav Hajek <highegg@gmail.com>
parents: 9628
diff changeset
234 s = - __doglegm__ (hesr, grad, dg, delta);
9084
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 sn = norm (dg .* s);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
237 if (niter == 1)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
238 delta = min (delta, sn);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
239 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
240
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
241 fval1 = fcn (reshape (x + s, xsiz)) (:);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
242 nfev ++;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
243
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
244 if (fval1 < fval)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
245 ## Scaled actual reduction.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
246 actred = (fval - fval1) / (abs (fval1) + abs (fval));
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
247 else
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
248 actred = -1;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
249 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
250
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
251 w = hesr*s;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
252 ## Scaled predicted reduction, and ratio.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
253 t = 1/2 * sumsq (w) + grad'*s;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
254 if (t < 0)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
255 prered = -t/(abs (fval) + abs (fval + t));
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
256 ratio = actred / prered;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
257 else
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
258 prered = 0;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
259 ratio = 0;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
260 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
261
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
262 ## Update delta.
9207
25f50d2d76b3 improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9199
diff changeset
263 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
264 delta *= decfac;
25f50d2d76b3 improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9199
diff changeset
265 decfac ^= 1.4142;
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
266 if (delta <= 1e1*macheps*xn)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
267 ## Trust region became uselessly small.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
268 info = -3;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
269 break;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
270 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
271 else
9199
399884c9d4a1 import the step adaptation strategy from fsolve to fminunc
Jaroslav Hajek <highegg@gmail.com>
parents: 9084
diff changeset
272 lastratio = ratio;
9207
25f50d2d76b3 improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9199
diff changeset
273 decfac = 0.5;
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
274 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
275 delta = 1.4142*sn;
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
276 elseif (ratio >= 0.5)
9207
25f50d2d76b3 improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents: 9199
diff changeset
277 delta = max (delta, 1.4142*sn);
9084
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 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
280
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
281 if (ratio >= 1e-4)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
282 ## Successful iteration.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
283 x += s;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
284 xn = norm (dg .* x);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
285 fval = fval1;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
286 nsuciter ++;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
287 suc = true;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
288 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
289
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
290 niter ++;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
291
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
292 ## 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
293 if (! isempty (outfcn))
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
294 optimvalues.iter = niter;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
295 optimvalues.funccount = nfev;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
296 optimvalues.fval = fval;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
297 optimvalues.searchdirection = s;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
298 state = 'iter';
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
299 stop = outfcn (x, optimvalues, state);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
300 if (stop)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
301 info = -1;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
302 break;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
303 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
304 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
305
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
306 ## Tests for termination conditions. A mysterious place, anything
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
307 ## 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
308
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
309 ## 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
310 ## 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
311 ## 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
312 ## iterations, so we need scaling-independent tolerances wherever
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
313 ## possible.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
314
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
315 ## The following tests done only after successful step.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
316 if (ratio >= 1e-4)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
317 ## 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
318 ## but compare to scaled step, so nothing bad.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
319 if (sn <= tolx*xn)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
320 info = 2;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
321 ## Again a classic one.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
322 elseif (actred < tolf)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
323 info = 3;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
324 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
325 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
326
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
327 endwhile
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
328 endwhile
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
329
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
330 ## Restore original shapes.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
331 x = reshape (x, xsiz);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
332
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
333 output.iterations = niter;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
334 output.successful = nsuciter;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
335 output.funcCount = nfev;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
336
9627
5bcfa0b346e8 fix extra outputs in fminunc
Jaroslav Hajek <highegg@gmail.com>
parents: 9623
diff changeset
337 if (nargout > 5)
5bcfa0b346e8 fix extra outputs in fminunc
Jaroslav Hajek <highegg@gmail.com>
parents: 9623
diff changeset
338 hess = hesr'*hesr;
5bcfa0b346e8 fix extra outputs in fminunc
Jaroslav Hajek <highegg@gmail.com>
parents: 9623
diff changeset
339 endif
5bcfa0b346e8 fix extra outputs in fminunc
Jaroslav Hajek <highegg@gmail.com>
parents: 9623
diff changeset
340
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
341 endfunction
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
342
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
343 ## 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
344 ## bad results.
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
345 function [fx, gx] = guarded_eval (fun, x)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
346 if (nargout > 1)
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
347 [fx, gx] = fun (x);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
348 else
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
349 fx = fun (x);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
350 gx = [];
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
351 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
352
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
353 if (! (isreal (fx) && isreal (jx)))
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
354 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
355 elseif (complexeqn && ! (isnumeric (fx) && isnumeric(jx)))
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
356 error ("fminunc:notnum", "fminunc: non-numeric value encountered");
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
357 elseif (any (isnan (fx(:))))
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
358 error ("fminunc:isnan", "fminunc: NaN value encountered");
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
359 endif
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
360 endfunction
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
361
13305
63463570d9fe Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents: 13027
diff changeset
362 %!function f = __rosenb (x)
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
363 %! n = length (x);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
364 %! 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
365 %!endfunction
9084
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
366 %!test
13305
63463570d9fe Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents: 13027
diff changeset
367 %! [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
368 %! tol = 2e-5;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
369 %! assert (info > 0);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
370 %! assert (x, ones (1, 2), tol);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
371 %! assert (fval, 0, tol);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
372 %!test
13305
63463570d9fe Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents: 13027
diff changeset
373 %! [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
374 %! tol = 2e-5;
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
375 %! assert (info > 0);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
376 %! assert (x, ones (1, 4), tol);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
377 %! assert (fval, 0, tol);
b7210faa3ed0 implement fminunc using factored trust-region BFGS
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
378
9899
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
379 ## 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
380 ## 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
381 ## 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
382
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
383 ## TODO: error checks
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
384 ## 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
385
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
386 function x = __doglegm__ (r, g, d, delta)
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
387 ## Get Gauss-Newton direction.
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
388 b = r' \ g;
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
389 x = r \ b;
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
390 xn = norm (d .* x);
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
391 if (xn > delta)
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
392 ## GN is too big, get scaled gradient.
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
393 s = g ./ d;
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
394 sn = norm (s);
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
395 if (sn > 0)
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
396 ## Normalize and rescale.
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
397 s = (s / sn) ./ d;
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
398 ## Get the line minimizer in s direction.
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
399 tn = norm (r*s);
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
400 snm = (sn / tn) / tn;
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
401 if (snm < delta)
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10201
diff changeset
402 ## Get the dogleg path minimizer.
9899
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
403 bn = norm (b);
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
404 dxn = delta/xn; snmd = snm/delta;
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
405 t = (bn/sn) * (bn/xn) * snmd;
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
406 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
407 alpha = dxn*(1-snmd^2) / t;
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
408 else
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
409 alpha = 0;
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
410 endif
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
411 else
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
412 alpha = delta / xn;
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
413 snm = 0;
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
414 endif
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
415 ## Form the appropriate convex combination.
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
416 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
417 endif
9f25290a35e8 more private function and subfunction changes
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
418 endfunction