Mercurial > hg > octave-nkf
annotate scripts/optimization/fminbnd.m @ 10297:ed88ea036716
improve docs of fzero/fminbnd
author | Jaroslav Hajek <highegg@gmail.com> |
---|---|
date | Wed, 10 Feb 2010 16:15:30 +0100 |
parents | 035ac548a67e |
children | b4e5dcf023c9 |
rev | line source |
---|---|
10296 | 1 ## Copyright (C) 2008, 2009, 2010 VZLU Prague, a.s. |
2 ## | |
3 ## This file is part of Octave. | |
4 ## | |
5 ## Octave is free software; you can redistribute it and/or modify it | |
6 ## under the terms of the GNU General Public License as published by | |
7 ## the Free Software Foundation; either version 3 of the License, or (at | |
8 ## your option) any later version. | |
9 ## | |
10 ## Octave is distributed in the hope that it will be useful, but | |
11 ## WITHOUT ANY WARRANTY; without even the implied warranty of | |
12 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU | |
13 ## General Public License for more details. | |
14 ## | |
15 ## You should have received a copy of the GNU General Public License | |
16 ## along with Octave; see the file COPYING. If not, see | |
17 ## <http://www.gnu.org/licenses/>. | |
18 ## | |
19 ## Author: Jaroslav Hajek <highegg@gmail.com> | |
20 | |
21 ## -*- texinfo -*- | |
22 ## @deftypefn {Function File} {[@var{x}, @var{fval}, @var{info}, @var{output}] =} fminbnd (@var{fun}, @var{a}, @var{b}, @var{options}) | |
23 ## Find a minimum point of a univariate function. @var{fun} should be a function | |
24 ## handle or name. @var{a}, @var{b} specify a starting interval. @var{options} is a | |
25 ## structure specifying additional options. Currently, @code{fminbnd} | |
26 ## recognizes these options: @code{"FunValCheck"}, @code{"OutputFcn"}, | |
27 ## @code{"TolX"}, @code{"MaxIter"}, @code{"MaxFunEvals"}. | |
28 ## For description of these options, see @ref{doc-optimset,,optimset}. | |
29 ## | |
30 ## On exit, the function returns @var{x}, the approximate minimum point | |
31 ## and @var{fval}, the function value thereof. | |
32 ## @var{info} is an exit flag that can have these values: | |
10297
ed88ea036716
improve docs of fzero/fminbnd
Jaroslav Hajek <highegg@gmail.com>
parents:
10296
diff
changeset
|
33 ## |
10296 | 34 ## @itemize |
35 ## @item 1 | |
36 ## The algorithm converged to a solution. | |
37 ## @item 0 | |
38 ## Maximum number of iterations or function evaluations has been exhausted. | |
39 ## @item -1 | |
40 ## The algorithm has been terminated from user output function. | |
41 ## @end itemize | |
42 ## @seealso{optimset, fzero, fminunc} | |
43 ## @end deftypefn | |
44 | |
45 ## This is patterned after opt/fmin.f from Netlib, which in turn is taken from | |
46 ## Richard Brent: Algorithms For Minimization Without Derivatives, Prentice-Hall (1973) | |
47 | |
48 ## PKG_ADD: __all_opts__ ("fminbnd"); | |
49 | |
50 function [x, fval, info, output] = fminbnd (fun, xmin, xmax, options = struct ()) | |
51 | |
52 ## Get default options if requested. | |
53 if (nargin == 1 && ischar (fun) && strcmp (fun, 'defaults')) | |
54 x = optimset ("MaxIter", Inf, "MaxFunEvals", Inf, "TolX", 1e-8, \ | |
55 "OutputFcn", [], "FunValCheck", "off"); | |
56 return; | |
57 endif | |
58 | |
59 if (nargin < 2 || nargin > 4) | |
60 print_usage (); | |
61 endif | |
62 | |
63 if (ischar (fun)) | |
64 fun = str2func (fun, "global"); | |
65 endif | |
66 | |
67 ## TODO | |
68 ## displev = optimget (options, "Display", "notify"); | |
69 funvalchk = strcmpi (optimget (options, "FunValCheck", "off"), "on"); | |
70 outfcn = optimget (options, "OutputFcn"); | |
71 tolx = optimget (options, "TolX", 1e-8); | |
72 maxiter = optimget (options, "MaxIter", Inf); | |
73 maxfev = optimget (options, "MaxFunEvals", Inf); | |
74 | |
75 persistent mu = 0.5; | |
76 | |
77 if (funvalchk) | |
78 ## Replace fun with a guarded version. | |
79 fun = @(x) guarded_eval (fun, x); | |
80 endif | |
81 | |
82 ## The default exit flag if exceeded number of iterations. | |
83 info = 0; | |
84 niter = 0; | |
85 nfev = 0; | |
86 eps = eps (class (xmin + xmax)); | |
87 | |
88 c = 0.5*(3-sqrt(5)); | |
89 a = xmin; b = xmax; | |
90 v = a + c*(b-a); | |
91 w = x = v; | |
92 e = 0; | |
93 fv = fw = fval = fun (x); | |
94 nfev++; | |
95 | |
96 while (niter < maxiter && nfev < maxfev) | |
97 xm = 0.5*(a+b); | |
98 tol = 2 * eps * abs (x) + tolx / 3; | |
99 if (abs (x - xm) <= (2*tol - 0.5*(b-a))) | |
100 info = 1; | |
101 break; | |
102 endif | |
103 | |
104 if (abs (e) > tol) | |
105 dogs = false; | |
106 ## Try inverse parabolic step. | |
107 r = (x - w)*(fval - fv); | |
108 q = (x - v)*(fval - fw); | |
109 p = (x - v)*q - (x - w)*r; | |
110 q = 2*(q - r); | |
111 p *= -sign (q); | |
112 q = abs (q); | |
113 r = e; | |
114 e = d; | |
115 | |
116 if (abs (p) < abs (0.5*q*r) && p > q*(a-x) && p < q*(b-x)) | |
117 ## The parabolic step is acceptable. | |
118 d = p / q; | |
119 u = x + d; | |
120 | |
121 ## f must not be evaluated too close to ax or bx. | |
122 if ((u-a) < 2*tol && (b-u) < 2*tol) | |
123 d = tol * (sign (xm - x) + (xm == x)); | |
124 endif | |
125 else | |
126 dogs = true; | |
127 endif | |
128 else | |
129 dogs = true; | |
130 endif | |
131 if (dogs) | |
132 ## Default to golden section step. | |
133 e = ifelse (x >= xm, a - x, b - x); | |
134 d = c * e; | |
135 endif | |
136 | |
137 ## f must not be evaluated too close to x. | |
138 u = x + max (abs (d), tol) * (sign (d) + (d == 0)); | |
139 | |
140 fu = fun (u); | |
141 nfev++; | |
142 niter++; | |
143 | |
144 ## update a, b, v, w, and x | |
145 | |
146 if (fu <= fval) | |
147 if (u < x) | |
148 b = x; | |
149 else | |
150 a = x; | |
151 endif | |
152 v = w; fv = fw; | |
153 w = x; fw = fval; | |
154 x = u; fval = fu; | |
155 else | |
156 ## The following if-statement was originally executed even if fu == fval. | |
157 if (u < x) | |
158 a = u; | |
159 else | |
160 b = u; | |
161 endif | |
162 if (fu <= fw || w == x) | |
163 v = w; fv = fw; | |
164 w = u; fw = fu; | |
165 elseif (fu <= fv || v == x || v == w) | |
166 v = u; | |
167 fv = fu; | |
168 endif | |
169 endif | |
170 | |
171 ## If there's an output function, use it now. | |
172 if (outfcn) | |
173 optv.funccount = nfev; | |
174 optv.fval = fval; | |
175 optv.iteration = niter; | |
176 if (outfcn (x, optv, "iter")) | |
177 info = -1; | |
178 break; | |
179 endif | |
180 endif | |
181 endwhile | |
182 | |
183 output.iterations = niter; | |
184 output.funcCount = nfev; | |
185 output.bracket = [a, b]; | |
186 ## FIXME: bracketf possibly unavailable. | |
187 | |
188 endfunction | |
189 | |
190 ## An assistant function that evaluates a function handle and checks for | |
191 ## bad results. | |
192 function fx = guarded_eval (fun, x) | |
193 fx = fun (x); | |
194 fx = fx(1); | |
195 if (! isreal (fx)) | |
196 error ("fminbnd:notreal", "fminbnd: non-real value encountered"); | |
197 elseif (isnan (fx)) | |
198 error ("fminbnd:isnan", "fminbnd: NaN value encountered"); | |
199 endif | |
200 endfunction | |
201 | |
202 %!shared opt0 | |
203 %! opt0 = optimset ("tolx", 0); | |
204 %!assert (fminbnd (@cos, pi/2, 3*pi/2, opt0), pi, 10*eps) | |
205 %!assert (fminbnd (@(x) (x - 1e-3)^4, -1, 1, opt0), 1e-3, 10e-3*eps) | |
206 %!assert (fminbnd (@(x) abs(x-1e7), 0, 1e10, opt0), 1e7, 10e7*eps) | |
207 %!assert (fminbnd (@(x) x^2 + sin(2*pi*x), 0.4, 1, opt0), fzero (@(x) 2*x + 2*pi*cos(2*pi*x), [0.4, 1], opt0), 1e-7) |