Mercurial > octave
annotate scripts/optimization/fsolve.m @ 33632:fed0dc6fd44c default tip
remove unused variable from libgui/module.mk
* libgui/module.mk: remove empty variable OCTAVE_GUI_EDITOR_MOC
author | Torsten Lilge <ttl-octave@mailbox.org> |
---|---|
date | Mon, 27 May 2024 19:42:05 +0200 |
parents | 2e484f9f1f18 |
children |
rev | line source |
---|---|
27923
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
1 ######################################################################## |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
2 ## |
32632
2e484f9f1f18
maint: update Octave Project Developers copyright for the new year
John W. Eaton <jwe@octave.org>
parents:
31706
diff
changeset
|
3 ## Copyright (C) 2008-2024 The Octave Project Developers |
27918
b442ec6dda5c
use centralized file for copyright info for individual contributors
John W. Eaton <jwe@octave.org>
parents:
27390
diff
changeset
|
4 ## |
27923
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
5 ## See the file COPYRIGHT.md in the top-level directory of this |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
6 ## distribution or <https://octave.org/copyright/>. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
7 ## |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
8 ## This file is part of Octave. |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
9 ## |
24534
194eb4bd202b
maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents:
23563
diff
changeset
|
10 ## Octave is free software: you can redistribute it and/or modify it |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
11 ## under the terms of the GNU General Public License as published by |
24534
194eb4bd202b
maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents:
23563
diff
changeset
|
12 ## the Free Software Foundation, either version 3 of the License, or |
22755
3a2b891d0b33
maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents:
22323
diff
changeset
|
13 ## (at your option) any later version. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
14 ## |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
15 ## Octave is distributed in the hope that it will be useful, but |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
16 ## WITHOUT ANY WARRANTY; without even the implied warranty of |
22755
3a2b891d0b33
maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents:
22323
diff
changeset
|
17 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
3a2b891d0b33
maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents:
22323
diff
changeset
|
18 ## GNU General Public License for more details. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
19 ## |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
20 ## You should have received a copy of the GNU General Public License |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
21 ## along with Octave; see the file COPYING. If not, see |
24534
194eb4bd202b
maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents:
23563
diff
changeset
|
22 ## <https://www.gnu.org/licenses/>. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
23 ## |
27923
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
24 ######################################################################## |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
25 |
8466 | 26 ## -*- texinfo -*- |
30875
5d3faba0342e
doc: Ensure documentation lists output argument when it exists for all m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
27 ## @deftypefn {} {@var{x} =} fsolve (@var{fcn}, @var{x0}) |
5d3faba0342e
doc: Ensure documentation lists output argument when it exists for all m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
28 ## @deftypefnx {} {@var{x} =} fsolve (@var{fcn}, @var{x0}, @var{options}) |
5d3faba0342e
doc: Ensure documentation lists output argument when it exists for all m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
29 ## @deftypefnx {} {[@var{x}, @var{fval}] =} fsolve (@dots{}) |
5d3faba0342e
doc: Ensure documentation lists output argument when it exists for all m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
30 ## @deftypefnx {} {[@var{x}, @var{fval}, @var{info}] =} fsolve (@dots{}) |
5d3faba0342e
doc: Ensure documentation lists output argument when it exists for all m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
31 ## @deftypefnx {} {[@var{x}, @var{fval}, @var{info}, @var{output}] =} fsolve (@dots{}) |
26075
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
32 ## @deftypefnx {} {[@var{x}, @var{fval}, @var{info}, @var{output}, @var{fjac}] =} fsolve (@dots{}) |
8513 | 33 ## Solve a system of nonlinear equations defined by the function @var{fcn}. |
20165
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
34 ## |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30875
diff
changeset
|
35 ## @var{fcn} is a function handle, inline function, or string containing the |
27390
baeed03c3766
doc: Use common verbiage to describe input FUN in optimization functions.
Rik <rik@octave.org>
parents:
27069
diff
changeset
|
36 ## name of the function to evaluate. @var{fcn} should accept a vector (array) |
baeed03c3766
doc: Use common verbiage to describe input FUN in optimization functions.
Rik <rik@octave.org>
parents:
27069
diff
changeset
|
37 ## defining the unknown variables, and return a vector of left-hand sides of |
baeed03c3766
doc: Use common verbiage to describe input FUN in optimization functions.
Rik <rik@octave.org>
parents:
27069
diff
changeset
|
38 ## the equations. Right-hand sides are defined to be zeros. In other words, |
baeed03c3766
doc: Use common verbiage to describe input FUN in optimization functions.
Rik <rik@octave.org>
parents:
27069
diff
changeset
|
39 ## this function attempts to determine a vector @var{x} such that |
baeed03c3766
doc: Use common verbiage to describe input FUN in optimization functions.
Rik <rik@octave.org>
parents:
27069
diff
changeset
|
40 ## @code{@var{fcn} (@var{x})} gives (approximately) all zeros. |
20165
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
41 ## |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
42 ## @var{x0} is an initial guess for the solution. The shape of @var{x0} is |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
43 ## preserved in all calls to @var{fcn}, but otherwise is treated as a column |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
44 ## vector. |
20165
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
45 ## |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
46 ## @var{options} is a structure specifying additional parameters which |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
47 ## control the algorithm. Currently, @code{fsolve} recognizes these options: |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
48 ## @qcode{"AutoScaling"}, @qcode{"ComplexEqn"}, @qcode{"FinDiffType"}, |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
49 ## @qcode{"FunValCheck"}, @qcode{"Jacobian"}, @qcode{"MaxFunEvals"}, |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
50 ## @qcode{"MaxIter"}, @qcode{"OutputFcn"}, @qcode{"TolFun"}, @qcode{"TolX"}, |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
51 ## @qcode{"TypicalX"}, and @qcode{"Updating"}. |
10201
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
52 ## |
26075
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
53 ## If @qcode{"AutoScaling"} is @qcode{"on"}, the variables will be |
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
54 ## automatically scaled according to the column norms of the (estimated) |
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
55 ## Jacobian. As a result, @qcode{"TolFun"} becomes scaling-independent. By |
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
56 ## default, this option is @qcode{"off"} because it may sometimes deliver |
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
57 ## unexpected (though mathematically correct) results. |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
58 ## |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
59 ## If @qcode{"ComplexEqn"} is @qcode{"on"}, @code{fsolve} will attempt to solve |
20165
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
60 ## complex equations in complex variables, assuming that the equations possess |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
61 ## a complex derivative (i.e., are holomorphic). If this is not what you want, |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
62 ## you should unpack the real and imaginary parts of the system to get a real |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
63 ## system. |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
64 ## |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
65 ## If @qcode{"Jacobian"} is @qcode{"on"}, it specifies that @var{fcn}---when |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
66 ## called with 2 output arguments---also returns the Jacobian matrix of |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
67 ## right-hand sides at the requested point. |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
68 ## |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
69 ## @qcode{"MaxFunEvals"} proscribes the maximum number of function evaluations |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
70 ## before optimization is halted. The default value is |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
71 ## @code{100 * number_of_variables}, i.e., @code{100 * length (@var{x0})}. |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
72 ## The value must be a positive integer. |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
73 ## |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
74 ## If @qcode{"Updating"} is @qcode{"on"}, the function will attempt to use |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
75 ## @nospell{Broyden} updates to update the Jacobian, in order to reduce the |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
76 ## number of Jacobian calculations. If your user function always calculates |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
77 ## the Jacobian (regardless of number of output arguments) then this option |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
78 ## provides no advantage and should be disabled. |
8466 | 79 ## |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
80 ## @qcode{"TolX"} specifies the termination tolerance in the unknown variables, |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
81 ## while @qcode{"TolFun"} is a tolerance for equations. Default is @code{1e-6} |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
82 ## for both @qcode{"TolX"} and @qcode{"TolFun"}. |
20165
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
83 ## |
28958
6fd6ad758b10
doc: Use @xref, @pxref rather than "see @code{}" in TexInfo.
Rik <rik@octave.org>
parents:
28789
diff
changeset
|
84 ## For a description of the other options, |
6fd6ad758b10
doc: Use @xref, @pxref rather than "see @code{}" in TexInfo.
Rik <rik@octave.org>
parents:
28789
diff
changeset
|
85 ## @pxref{XREFoptimset,,@code{optimset}}. To initialize an options structure |
6fd6ad758b10
doc: Use @xref, @pxref rather than "see @code{}" in TexInfo.
Rik <rik@octave.org>
parents:
28789
diff
changeset
|
86 ## with default values for @code{fsolve} use |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
87 ## @code{options = optimset ("fsolve")}. |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
88 ## |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
89 ## The first output @var{x} is the solution while the second output @var{fval} |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
90 ## contains the value of the function @var{fcn} evaluated at @var{x} (ideally |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
91 ## a vector of all zeros). |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
92 ## |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
93 ## The third output @var{info} reports whether the algorithm succeeded and may |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
94 ## take 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
|
95 ## |
8466 | 96 ## @table @asis |
97 ## @item 1 | |
10793
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
10711
diff
changeset
|
98 ## Converged to a solution point. Relative residual error is less than |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
99 ## specified by @code{TolFun}. |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
100 ## |
8466 | 101 ## @item 2 |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
102 ## Last relative step size was less than @code{TolX}. |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
103 ## |
8466 | 104 ## @item 3 |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
105 ## Last relative decrease in residual was less than @code{TolFun}. |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
106 ## |
8466 | 107 ## @item 0 |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
108 ## Iteration limit (either @code{MaxIter} or @code{MaxFunEvals}) exceeded. |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
109 ## |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
110 ## @item -1 |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
111 ## Stopped by @code{OutputFcn}. |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
112 ## |
30234
3dfe97a06ba9
doc: Document fsolve output "info" -2 (bug #61310).
Arun Giridhar <arungiridhar@gmail.com>
parents:
29397
diff
changeset
|
113 ## @item -2 |
3dfe97a06ba9
doc: Document fsolve output "info" -2 (bug #61310).
Arun Giridhar <arungiridhar@gmail.com>
parents:
29397
diff
changeset
|
114 ## The Jacobian became excessively small and the search stalled. |
3dfe97a06ba9
doc: Document fsolve output "info" -2 (bug #61310).
Arun Giridhar <arungiridhar@gmail.com>
parents:
29397
diff
changeset
|
115 ## |
8466 | 116 ## @item -3 |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
117 ## The trust region radius became excessively small. |
8466 | 118 ## @end table |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
119 ## |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
120 ## @var{output} is a structure containing runtime information about the |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
121 ## @code{fsolve} algorithm. Fields in the structure are: |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
122 ## |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
123 ## @table @code |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
124 ## @item iterations |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
125 ## Number of iterations through loop. |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
126 ## |
26075
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
127 ## @item successful |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
128 ## Number of successful iterations. |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
129 ## |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
130 ## @item @nospell{funcCount} |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
131 ## Number of function evaluations. |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
132 ## |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
133 ## @end table |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
134 ## |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
135 ## The final output @var{fjac} contains the value of the Jacobian evaluated |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
136 ## at @var{x}. |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
137 ## |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
138 ## Note: If you only have a single nonlinear equation of one variable, using |
8466 | 139 ## @code{fzero} is usually a much better idea. |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
140 ## |
10793
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
10711
diff
changeset
|
141 ## Note about user-supplied Jacobians: |
20165
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
142 ## As an inherent property of the algorithm, a Jacobian is always requested for |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
143 ## a solution vector whose residual vector is already known, and it is the last |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8997
diff
changeset
|
144 ## accepted successful step. Often this will be one of the last two calls, but |
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8997
diff
changeset
|
145 ## not always. If the savings by reusing intermediate results from residual |
10793
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
10711
diff
changeset
|
146 ## calculation in Jacobian calculation are significant, the best strategy is to |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
147 ## employ @code{OutputFcn}: After a vector is evaluated for residuals, if |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
148 ## @code{OutputFcn} is called with that vector, then the intermediate results |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
149 ## should be saved for future Jacobian evaluation, and should be kept until a |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
150 ## Jacobian evaluation is requested or until @code{OutputFcn} is called with a |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
151 ## different vector, in which case they should be dropped in favor of this most |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
152 ## recent vector. A short example how this can be achieved follows: |
8988
315828058e0d
minor doc fixes in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8986
diff
changeset
|
153 ## |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
154 ## @example |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30875
diff
changeset
|
155 ## function [fval, fjac] = user_fcn (x, optimvalues, state) |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
156 ## persistent sav = [], sav0 = []; |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
157 ## if (nargin == 1) |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
158 ## ## evaluation call |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
159 ## if (nargout == 1) |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
160 ## sav0.x = x; # mark saved vector |
26075
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
161 ## ## calculate fval, save results to sav0. |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
162 ## elseif (nargout == 2) |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
163 ## ## calculate fjac using sav. |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
164 ## endif |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
165 ## else |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
166 ## ## outputfcn call. |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
167 ## if (all (x == sav0.x)) |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
168 ## sav = sav0; |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
169 ## endif |
9209
923c7cb7f13f
Simplify TeXinfo files by eliminating redundant @iftex followed by @tex construction.
Rik <rdrider0-list@yahoo.com>
parents:
9207
diff
changeset
|
170 ## ## maybe output iteration status, etc. |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
171 ## endif |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
172 ## endfunction |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
173 ## |
14143
04dcbb8fb880
fsolve.m: Move @seealso to bottom in docstring to silence warning.
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
14138
diff
changeset
|
174 ## ## @dots{} |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
175 ## |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30875
diff
changeset
|
176 ## fsolve (@@user_fcn, x0, optimset ("OutputFcn", @@user_fcn, @dots{})) |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
177 ## @end example |
14143
04dcbb8fb880
fsolve.m: Move @seealso to bottom in docstring to silence warning.
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
14138
diff
changeset
|
178 ## @seealso{fzero, optimset} |
8466 | 179 ## @end deftypefn |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
180 |
13027
b9a89ca0fb75
prevent optimization functions from setting ans in workspace at startup
John W. Eaton <jwe@octave.org>
parents:
11587
diff
changeset
|
181 ## 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:
11587
diff
changeset
|
182 ## PKG_ADD: [~] = __all_opts__ ("fsolve"); |
8648
ff61b53eb294
optimization: use PKG_ADD: comments instead of PKG_ADD file
John W. Eaton <jwe@octave.org>
parents:
8647
diff
changeset
|
183 |
26075
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
184 function [x, fval, info, output, fjac] = fsolve (fcn, x0, options = struct ()) |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
185 |
8596
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
186 ## Get default options if requested. |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
187 if (nargin == 1 && ischar (fcn) && strcmp (fcn, "defaults")) |
26138
804e18e3e320
Reenable query of optimization options (bugs #54952 and #55089).
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
26118
diff
changeset
|
188 x = struct ("AutoScaling", "off", "ComplexEqn", "off", |
804e18e3e320
Reenable query of optimization options (bugs #54952 and #55089).
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
26118
diff
changeset
|
189 "FunValCheck", "off", "FinDiffType", "forward", |
804e18e3e320
Reenable query of optimization options (bugs #54952 and #55089).
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
26118
diff
changeset
|
190 "Jacobian", "off", "MaxFunEvals", [], "MaxIter", 400, |
804e18e3e320
Reenable query of optimization options (bugs #54952 and #55089).
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
26118
diff
changeset
|
191 "OutputFcn", [], "Updating", "off", "TolFun", 1e-6, |
804e18e3e320
Reenable query of optimization options (bugs #54952 and #55089).
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
26118
diff
changeset
|
192 "TolX", 1e-6, "TypicalX", []); |
8596
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
193 return; |
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
194 endif |
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
195 |
28789
28de41192f3c
Eliminate unneeded verification of nargin, nargout in m-files.
Rik <rik@octave.org>
parents:
27978
diff
changeset
|
196 if (nargin < 2 || ! isnumeric (x0)) |
8596
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
197 print_usage (); |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
198 endif |
8596
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
199 |
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
200 if (ischar (fcn)) |
26118
7502fce4cd3a
str2func: eliminate optional second "global" argument
John W. Eaton <jwe@octave.org>
parents:
26075
diff
changeset
|
201 fcn = str2func (fcn); |
10050
dc88a0b6472c
support old style jacobian for fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
202 elseif (iscell (fcn)) |
dc88a0b6472c
support old style jacobian for fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
203 fcn = @(x) make_fcn_jac (x, fcn{1}, fcn{2}); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
204 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
205 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
206 xsiz = size (x0); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
207 n = numel (x0); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
208 |
8467
77b8d4aa2743
fsolve.m, fzero.m: style fixes; use strcmpi to compare options
John W. Eaton <jwe@octave.org>
parents:
8466
diff
changeset
|
209 has_jac = strcmpi (optimget (options, "Jacobian", "off"), "on"); |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
210 cdif = strcmpi (optimget (options, "FinDiffType", "forward"), "central"); |
8596
8833c0b18eb2
enable default settings queries in optim funcs
Jaroslav Hajek <highegg@gmail.com>
parents:
8592
diff
changeset
|
211 maxiter = optimget (options, "MaxIter", 400); |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
212 maxfev = optimget (options, "MaxFunEvals", 100*n); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
213 outfcn = optimget (options, "OutputFcn"); |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
214 updating = strcmpi (optimget (options, "Updating", "off"), "on"); |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
215 complexeqn = strcmpi (optimget (options, "ComplexEqn", "off"), "on"); |
8590 | 216 |
21751
b571fc85953f
maint: Use two spaces after period to indicate sentence break.
Rik <rik@octave.org>
parents:
21578
diff
changeset
|
217 ## Get scaling matrix using the TypicalX option. If set to "auto", the |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
218 ## scaling matrix is estimated using the Jacobian. |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
219 typicalx = optimget (options, "TypicalX", ones (n, 1)); |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
220 |
10201
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
221 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:
10050
diff
changeset
|
222 if (! autoscale) |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
223 dg = 1 ./ typicalx; |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
224 endif |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
225 |
8467
77b8d4aa2743
fsolve.m, fzero.m: style fixes; use strcmpi to compare options
John W. Eaton <jwe@octave.org>
parents:
8466
diff
changeset
|
226 funvalchk = strcmpi (optimget (options, "FunValCheck", "off"), "on"); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
227 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
228 if (funvalchk) |
8604 | 229 ## Replace fcn with a guarded version. |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
230 fcn = @(x) guarded_eval (fcn, x, complexeqn); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
231 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
232 |
21099
52af4092f863
For optimization scripts, correctly choose tolerance (eps) based on class of fun and X0.
Rik <rik@octave.org>
parents:
20852
diff
changeset
|
233 ## These defaults are rather stringent. |
52af4092f863
For optimization scripts, correctly choose tolerance (eps) based on class of fun and X0.
Rik <rik@octave.org>
parents:
20852
diff
changeset
|
234 ## Normally user prefers accuracy to performance. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
235 |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
236 tolx = optimget (options, "TolX", 1e-6); |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
237 tolf = optimget (options, "TolFun", 1e-6); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
238 |
9623
bc0739d02724
update initial TR step for fsolve and fminunc
Jaroslav Hajek <highegg@gmail.com>
parents:
9464
diff
changeset
|
239 factor = 1; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
240 |
8467
77b8d4aa2743
fsolve.m, fzero.m: style fixes; use strcmpi to compare options
John W. Eaton <jwe@octave.org>
parents:
8466
diff
changeset
|
241 niter = 1; |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
242 nfev = 1; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
243 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
244 x = x0(:); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
245 info = 0; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
246 |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
247 ## Initial evaluation. |
8990
349d75161672
more cosmetic adjustments to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8988
diff
changeset
|
248 ## Handle arbitrary shapes of x and f and remember them. |
26075
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
249 fval = fcn (reshape (x, xsiz)); |
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
250 fsiz = size (fval); |
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
251 fval = fval(:); |
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
252 fn = norm (fval); |
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
253 m = length (fval); |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
254 n = length (x); |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
255 |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
256 if (! isempty (outfcn)) |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
257 optimvalues.iter = niter; |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
258 optimvalues.funccount = nfev; |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
259 optimvalues.fval = fn; |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
260 optimvalues.searchdirection = zeros (n, 1); |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
261 state = "init"; |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
262 stop = outfcn (x, optimvalues, state); |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
263 if (stop) |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
264 info = -1; |
29397
ce425b657693
fsolve.m: Fix undefined output error when using Output function (bug #60144).
Jim Van Zandt <jim.vanzandt@gmail.com>
parents:
29358
diff
changeset
|
265 output.iterations = niter; |
ce425b657693
fsolve.m: Fix undefined output error when using Output function (bug #60144).
Jim Van Zandt <jim.vanzandt@gmail.com>
parents:
29358
diff
changeset
|
266 output.successful = 0; |
ce425b657693
fsolve.m: Fix undefined output error when using Output function (bug #60144).
Jim Van Zandt <jim.vanzandt@gmail.com>
parents:
29358
diff
changeset
|
267 output.funcCount = nfev; |
ce425b657693
fsolve.m: Fix undefined output error when using Output function (bug #60144).
Jim Van Zandt <jim.vanzandt@gmail.com>
parents:
29358
diff
changeset
|
268 fjac = NaN; |
22785
9c6661004167
error if break statement is in script file separate from loop (bug #39168)
John W. Eaton <jwe@octave.org>
parents:
22755
diff
changeset
|
269 return; |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
270 endif |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
271 endif |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
272 |
26075
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
273 if (isa (x0, "single") || isa (fval, "single")) |
21099
52af4092f863
For optimization scripts, correctly choose tolerance (eps) based on class of fun and X0.
Rik <rik@octave.org>
parents:
20852
diff
changeset
|
274 macheps = eps ("single"); |
52af4092f863
For optimization scripts, correctly choose tolerance (eps) based on class of fun and X0.
Rik <rik@octave.org>
parents:
20852
diff
changeset
|
275 else |
52af4092f863
For optimization scripts, correctly choose tolerance (eps) based on class of fun and X0.
Rik <rik@octave.org>
parents:
20852
diff
changeset
|
276 macheps = eps ("double"); |
52af4092f863
For optimization scripts, correctly choose tolerance (eps) based on class of fun and X0.
Rik <rik@octave.org>
parents:
20852
diff
changeset
|
277 endif |
52af4092f863
For optimization scripts, correctly choose tolerance (eps) based on class of fun and X0.
Rik <rik@octave.org>
parents:
20852
diff
changeset
|
278 |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
279 nsuciter = 0; |
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
280 |
8507 | 281 ## Outer loop. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
282 while (niter < maxiter && nfev < maxfev && ! info) |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
283 |
8507 | 284 ## Calculate function value and Jacobian (possibly via FD). |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
285 if (has_jac) |
26075
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
286 [fval, fjac] = fcn (reshape (x, xsiz)); |
23563
614d71cdf614
maint: Strip trailing whitespace from files.
John W. Eaton <jwe@octave.org>
parents:
23474
diff
changeset
|
287 if (! all (size (fjac) == [m, n])) |
23474
70edb5512c6e
fsolve.m: Issue an error for incorrectly sized user Jacobian (bug #50454).
Rik <rik@octave.org>
parents:
23220
diff
changeset
|
288 error ("fsolve: Jacobian size should be (%d,%d), not (%d,%d)", |
70edb5512c6e
fsolve.m: Issue an error for incorrectly sized user Jacobian (bug #50454).
Rik <rik@octave.org>
parents:
23220
diff
changeset
|
289 m, n, rows (fjac), columns (fjac)); |
70edb5512c6e
fsolve.m: Issue an error for incorrectly sized user Jacobian (bug #50454).
Rik <rik@octave.org>
parents:
23220
diff
changeset
|
290 endif |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
291 ## If the Jacobian is sparse, disable Broyden updating. |
8592
dacfd030633a
handle sparse jacobians in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8590
diff
changeset
|
292 if (issparse (fjac)) |
dacfd030633a
handle sparse jacobians in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8590
diff
changeset
|
293 updating = false; |
dacfd030633a
handle sparse jacobians in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8590
diff
changeset
|
294 endif |
26075
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
295 fval = fval(:); |
20735
418ae0cb752f
Replace ++,-- with in-place operators for performance.
Rik <rik@octave.org>
parents:
20714
diff
changeset
|
296 nfev += 1; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
297 else |
26075
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
298 fjac = __fdjac__ (fcn, reshape (x, xsiz), fval, typicalx, cdif); |
9212
6feb27c38da1
support central differences in fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9209
diff
changeset
|
299 nfev += (1 + cdif) * length (x); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
300 endif |
8590 | 301 |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
302 ## For square and overdetermined systems, we update a QR factorization of |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
303 ## the Jacobian to avoid solving a full system in each step. In this case, |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
304 ## we pass a triangular matrix to __dogleg__. |
8590 | 305 useqr = updating && m >= n && n > 10; |
306 | |
307 if (useqr) | |
8616
3d75d717cbe0
do not pivot by default in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8615
diff
changeset
|
308 ## FIXME: Currently, pivoting is mostly useless because the \ operator |
3d75d717cbe0
do not pivot by default in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8615
diff
changeset
|
309 ## cannot exploit the resulting props of the triangular factor. |
3d75d717cbe0
do not pivot by default in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8615
diff
changeset
|
310 ## Unpivoted QR is significantly faster so it doesn't seem right to pivot |
21937
55f7de37b618
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
21759
diff
changeset
|
311 ## just to get invariance. Original MINPACK didn't pivot either, |
21751
b571fc85953f
maint: Use two spaces after period to indicate sentence break.
Rik <rik@octave.org>
parents:
21578
diff
changeset
|
312 ## at least when qr updating was used. |
8616
3d75d717cbe0
do not pivot by default in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8615
diff
changeset
|
313 [q, r] = qr (fjac, 0); |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
314 endif |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
315 |
10201
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
316 if (autoscale) |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
317 ## Get column norms, use them as scaling factors. |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
318 jcn = norm (fjac, "columns").'; |
10201
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
319 if (niter == 1) |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
320 dg = jcn; |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
321 dg(dg == 0) = 1; |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
322 else |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
323 ## Rescale adaptively. |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
324 ## FIXME: the original minpack used the following rescaling strategy: |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
325 ## dg = max (dg, jcn); |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
326 ## but it seems not good if we start with a bad guess yielding Jacobian |
19833
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
327 ## columns with large norms that later decrease, because the |
21751
b571fc85953f
maint: Use two spaces after period to indicate sentence break.
Rik <rik@octave.org>
parents:
21578
diff
changeset
|
328 ## corresponding variable will still be overscaled. Instead, we only |
19833
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
329 ## give the old scaling a small momentum, but do not honor it. |
10201
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
330 |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
331 dg = max (0.1*dg, jcn); |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
332 endif |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
333 endif |
5c66978f3fdf
support TypicalX and AutoScaling in fsolve/fminunc, don't autoscale by default
Jaroslav Hajek <highegg@gmail.com>
parents:
10050
diff
changeset
|
334 |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
335 if (niter == 1) |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
336 xn = norm (dg .* x); |
8590 | 337 ## FIXME: something better? |
9628
73e6ad869f08
further correct initial TR step strategy
Jaroslav Hajek <highegg@gmail.com>
parents:
9623
diff
changeset
|
338 delta = factor * max (xn, 1); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
339 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
340 |
8990
349d75161672
more cosmetic adjustments to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8988
diff
changeset
|
341 ## It also seems that in the case of fast (and inhomogeneously) changing |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
342 ## Jacobian, the Broyden updates are of little use, so maybe we could |
21751
b571fc85953f
maint: Use two spaces after period to indicate sentence break.
Rik <rik@octave.org>
parents:
21578
diff
changeset
|
343 ## skip them if a big disproportional change is expected. The question is, |
8990
349d75161672
more cosmetic adjustments to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8988
diff
changeset
|
344 ## of course, how to define the above terms :) |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
345 |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
346 lastratio = 0; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
347 nfail = 0; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
348 nsuc = 0; |
9207
25f50d2d76b3
improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9153
diff
changeset
|
349 decfac = 0.5; |
25f50d2d76b3
improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9153
diff
changeset
|
350 |
8507 | 351 ## Inner loop. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
352 while (niter <= maxiter && nfev < maxfev && ! info) |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
353 |
8590 | 354 ## Get trust-region model (dogleg) minimizer. |
355 if (useqr) | |
25486
dffd9f6ee85c
fsolve.m: exit with info=-2 when singularity reached (bug #53991).
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
356 if (norm (r, 1) < macheps * rows (r)) |
dffd9f6ee85c
fsolve.m: exit with info=-2 when singularity reached (bug #53991).
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
357 info = -2; |
dffd9f6ee85c
fsolve.m: exit with info=-2 when singularity reached (bug #53991).
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
358 break; |
dffd9f6ee85c
fsolve.m: exit with info=-2 when singularity reached (bug #53991).
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
359 endif |
26075
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
360 qtf = q'*fval; |
8590 | 361 s = - __dogleg__ (r, qtf, dg, delta); |
362 w = qtf + r * s; | |
363 else | |
25486
dffd9f6ee85c
fsolve.m: exit with info=-2 when singularity reached (bug #53991).
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
364 if (norm (fjac, 1) < macheps * rows (fjac)) |
dffd9f6ee85c
fsolve.m: exit with info=-2 when singularity reached (bug #53991).
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
365 info = -2; |
dffd9f6ee85c
fsolve.m: exit with info=-2 when singularity reached (bug #53991).
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
366 break; |
dffd9f6ee85c
fsolve.m: exit with info=-2 when singularity reached (bug #53991).
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
367 endif |
26075
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
368 s = - __dogleg__ (fjac, fval, dg, delta); |
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
369 w = fval + fjac * s; |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
370 endif |
8590 | 371 |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
372 sn = norm (dg .* s); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
373 if (niter == 1) |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
374 delta = min (delta, sn); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
375 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
376 |
26075
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
377 fval1 = fcn (reshape (x + s, xsiz)) (:); |
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
378 fn1 = norm (fval1); |
20735
418ae0cb752f
Replace ++,-- with in-place operators for performance.
Rik <rik@octave.org>
parents:
20714
diff
changeset
|
379 nfev += 1; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
380 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
381 if (fn1 < fn) |
8507 | 382 ## Scaled actual reduction. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
383 actred = 1 - (fn1/fn)^2; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
384 else |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
385 actred = -1; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
386 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
387 |
8507 | 388 ## Scaled predicted reduction, and ratio. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
389 t = norm (w); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
390 if (t < fn) |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
391 prered = 1 - (t/fn)^2; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
392 ratio = actred / prered; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
393 else |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
394 prered = 0; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
395 ratio = 0; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
396 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
397 |
8507 | 398 ## Update delta. |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14386
diff
changeset
|
399 if (ratio < min (max (0.1, 0.8*lastratio), 0.9)) |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
400 nsuc = 0; |
20735
418ae0cb752f
Replace ++,-- with in-place operators for performance.
Rik <rik@octave.org>
parents:
20714
diff
changeset
|
401 nfail += 1; |
9207
25f50d2d76b3
improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9153
diff
changeset
|
402 delta *= decfac; |
25f50d2d76b3
improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9153
diff
changeset
|
403 decfac ^= 1.4142; |
25500
b3c35a130f94
fsolve.m: Return info=1 when initial guess (x0) is correct (bug #53991).
Rik <rik@octave.org>
parents:
25499
diff
changeset
|
404 if (fn <= tolf*n*xn) |
b3c35a130f94
fsolve.m: Return info=1 when initial guess (x0) is correct (bug #53991).
Rik <rik@octave.org>
parents:
25499
diff
changeset
|
405 info = 1; |
b3c35a130f94
fsolve.m: Return info=1 when initial guess (x0) is correct (bug #53991).
Rik <rik@octave.org>
parents:
25499
diff
changeset
|
406 elseif (delta <= 1e1*macheps*xn) |
8507 | 407 ## Trust region became uselessly small. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
408 info = -3; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
409 break; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
410 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
411 else |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
412 lastratio = ratio; |
9207
25f50d2d76b3
improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9153
diff
changeset
|
413 decfac = 0.5; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
414 nfail = 0; |
20735
418ae0cb752f
Replace ++,-- with in-place operators for performance.
Rik <rik@octave.org>
parents:
20714
diff
changeset
|
415 nsuc += 1; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
416 if (abs (1-ratio) <= 0.1) |
9207
25f50d2d76b3
improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9153
diff
changeset
|
417 delta = 1.4142*sn; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
418 elseif (ratio >= 0.5 || nsuc > 1) |
9207
25f50d2d76b3
improve TR updating strategy for fminunc and fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9153
diff
changeset
|
419 delta = max (delta, 1.4142*sn); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
420 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
421 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
422 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
423 if (ratio >= 1e-4) |
8507 | 424 ## Successful iteration. |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
425 x += s; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
426 xn = norm (dg .* x); |
26075
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
427 fval = fval1; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
428 fn = fn1; |
20735
418ae0cb752f
Replace ++,-- with in-place operators for performance.
Rik <rik@octave.org>
parents:
20714
diff
changeset
|
429 nsuciter += 1; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
430 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
431 |
20735
418ae0cb752f
Replace ++,-- with in-place operators for performance.
Rik <rik@octave.org>
parents:
20714
diff
changeset
|
432 niter += 1; |
8857 | 433 |
8590 | 434 ## FIXME: should outputfcn be only called after a successful iteration? |
8615 | 435 if (! isempty (outfcn)) |
8590 | 436 optimvalues.iter = niter; |
437 optimvalues.funccount = nfev; | |
438 optimvalues.fval = fn; | |
439 optimvalues.searchdirection = s; | |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
440 state = "iter"; |
8590 | 441 stop = outfcn (x, optimvalues, state); |
442 if (stop) | |
443 info = -1; | |
29397
ce425b657693
fsolve.m: Fix undefined output error when using Output function (bug #60144).
Jim Van Zandt <jim.vanzandt@gmail.com>
parents:
29358
diff
changeset
|
444 output.iterations = niter; |
ce425b657693
fsolve.m: Fix undefined output error when using Output function (bug #60144).
Jim Van Zandt <jim.vanzandt@gmail.com>
parents:
29358
diff
changeset
|
445 output.successful = nsuciter; |
ce425b657693
fsolve.m: Fix undefined output error when using Output function (bug #60144).
Jim Van Zandt <jim.vanzandt@gmail.com>
parents:
29358
diff
changeset
|
446 output.funcCount = nfev; |
8590 | 447 break; |
448 endif | |
449 endif | |
450 | |
21751
b571fc85953f
maint: Use two spaces after period to indicate sentence break.
Rik <rik@octave.org>
parents:
21578
diff
changeset
|
451 ## Tests for termination conditions. A mysterious place, anything |
8466 | 452 ## 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
|
453 |
8466 | 454 ## The rule of thumb (which I'm not sure M*b is quite following) |
455 ## is that for a tolerance that depends on scaling, only 0 makes | |
21751
b571fc85953f
maint: Use two spaces after period to indicate sentence break.
Rik <rik@octave.org>
parents:
21578
diff
changeset
|
456 ## sense as a default value. But 0 usually means uselessly long |
8466 | 457 ## iterations, so we need scaling-independent tolerances wherever |
458 ## possible. | |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
459 |
18857
7bbe3658c5ef
maint: Use "FIXME:" coding convention in m-files.
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
460 ## FIXME: Why tolf*n*xn? If abs (e) ~ abs(x) * eps is a vector |
21751
b571fc85953f
maint: Use two spaces after period to indicate sentence break.
Rik <rik@octave.org>
parents:
21578
diff
changeset
|
461 ## of perturbations of x, then norm (fjac*e) <= eps*n*xn, i.e., by |
8466 | 462 ## tolf ~ eps we demand as much accuracy as we can expect. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
463 if (fn <= tolf*n*xn) |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
464 info = 1; |
10549 | 465 ## The following tests done only after successful step. |
9075 | 466 elseif (ratio >= 1e-4) |
21751
b571fc85953f
maint: Use two spaces after period to indicate sentence break.
Rik <rik@octave.org>
parents:
21578
diff
changeset
|
467 ## This one is classic. Note that we use scaled variables again, |
10549 | 468 ## but compare to scaled step, so nothing bad. |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
469 if (sn <= tolx*xn) |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
470 info = 2; |
21751
b571fc85953f
maint: Use two spaces after period to indicate sentence break.
Rik <rik@octave.org>
parents:
21578
diff
changeset
|
471 ## Again a classic one. It seems weird to use the same tolf |
10549 | 472 ## for two different tests, but that's what M*b manual appears |
473 ## to say. | |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
474 elseif (actred < tolf) |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
475 info = 3; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
476 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
477 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
478 |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
479 ## Criterion for recalculating Jacobian. |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
480 if (! updating || nfail == 2 || nsuciter < 2) |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
481 break; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
482 endif |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
483 |
8507 | 484 ## Compute the scaled Broyden update. |
8590 | 485 if (useqr) |
26075
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
486 u = (fval1 - q*w) / sn; |
8590 | 487 v = dg .* ((dg .* s) / sn); |
488 | |
489 ## Update the QR factorization. | |
490 [q, r] = qrupdate (q, r, u, v); | |
491 else | |
26075
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
492 u = (fval1 - w); |
8590 | 493 v = dg .* ((dg .* s) / sn); |
494 | |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
495 ## update the Jacobian |
8590 | 496 fjac += u * v'; |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
497 endif |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
498 endwhile |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
499 endwhile |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
500 |
8507 | 501 ## Restore original shapes. |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
502 x = reshape (x, xsiz); |
26075
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
503 fval = reshape (fval, fsiz); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
504 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
505 output.iterations = niter; |
8986
22c8272af34b
improvements to fsolve & family
Jaroslav Hajek <highegg@gmail.com>
parents:
8857
diff
changeset
|
506 output.successful = nsuciter; |
8590 | 507 output.funcCount = nfev; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
508 |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
509 endfunction |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
510 |
21759
b002b4331a12
maint: Use two newlines after endfunction and start of BIST tests.
Rik <rik@octave.org>
parents:
21758
diff
changeset
|
511 ## A helper function that evaluates a function and checks for bad results. |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30875
diff
changeset
|
512 function [fx, jx] = guarded_eval (fcn, x, complexeqn) |
21758
ffad2baa90f7
maint: Use newlines to make code more readable.
Rik <rik@octave.org>
parents:
21751
diff
changeset
|
513 |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
514 if (nargout > 1) |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30875
diff
changeset
|
515 [fx, jx] = fcn (x); |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
516 else |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30875
diff
changeset
|
517 fx = fcn (x); |
8997
187a9d9c2f04
simplifications to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8991
diff
changeset
|
518 jx = []; |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
519 endif |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
520 |
8997
187a9d9c2f04
simplifications to fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8991
diff
changeset
|
521 if (! complexeqn && ! (isreal (fx) && isreal (jx))) |
27069
0a62d9a6aa2d
Place Octave's warning and error IDs in to the "Octave" namespace (bug #56213).
Rik <rik@octave.org>
parents:
26376
diff
changeset
|
522 error ("Octave:fsolve:notreal", "fsolve: non-real value encountered"); |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14386
diff
changeset
|
523 elseif (complexeqn && ! (isnumeric (fx) && isnumeric (jx))) |
27069
0a62d9a6aa2d
Place Octave's warning and error IDs in to the "Octave" namespace (bug #56213).
Rik <rik@octave.org>
parents:
26376
diff
changeset
|
524 error ("Octave:fsolve:notnum", "fsolve: non-numeric value encountered"); |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
525 elseif (any (isnan (fx(:)))) |
27069
0a62d9a6aa2d
Place Octave's warning and error IDs in to the "Octave" namespace (bug #56213).
Rik <rik@octave.org>
parents:
26376
diff
changeset
|
526 error ("Octave:fsolve:isnan", "fsolve: NaN value encountered"); |
14386
59aab666f2bf
Extend "FunValCheck" option to optimization routines to detect Inf values.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
527 elseif (any (isinf (fx(:)))) |
27069
0a62d9a6aa2d
Place Octave's warning and error IDs in to the "Octave" namespace (bug #56213).
Rik <rik@octave.org>
parents:
26376
diff
changeset
|
528 error ("Octave:fsolve:isinf", "fsolve: Inf value encountered"); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
529 endif |
21758
ffad2baa90f7
maint: Use newlines to make code more readable.
Rik <rik@octave.org>
parents:
21751
diff
changeset
|
530 |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
531 endfunction |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
532 |
10050
dc88a0b6472c
support old style jacobian for fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
533 function [fx, jx] = make_fcn_jac (x, fcn, fjac) |
21758
ffad2baa90f7
maint: Use newlines to make code more readable.
Rik <rik@octave.org>
parents:
21751
diff
changeset
|
534 |
10050
dc88a0b6472c
support old style jacobian for fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
535 fx = fcn (x); |
dc88a0b6472c
support old style jacobian for fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
536 if (nargout == 2) |
dc88a0b6472c
support old style jacobian for fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
537 jx = fjac (x); |
dc88a0b6472c
support old style jacobian for fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
538 endif |
21758
ffad2baa90f7
maint: Use newlines to make code more readable.
Rik <rik@octave.org>
parents:
21751
diff
changeset
|
539 |
10050
dc88a0b6472c
support old style jacobian for fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
540 endfunction |
dc88a0b6472c
support old style jacobian for fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
9899
diff
changeset
|
541 |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
542 ## Solve the double dogleg trust-region least-squares problem: |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
543 ## Minimize norm (r*x-b) subject to the constraint norm (d.*x) <= delta, |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
544 ## x being a convex combination of the gauss-newton and scaled gradient. |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
545 |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
546 ## FIXME: error checks |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
547 ## FIXME: handle singularity, or leave it up to mldivide? |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
548 |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
549 function x = __dogleg__ (r, b, d, delta) |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
550 |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
551 ## Get Gauss-Newton direction. |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
552 x = r \ b; |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
553 xn = norm (d .* x); |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
554 if (xn > delta) |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
555 ## GN is too big, get scaled gradient. |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
556 s = (r' * b) ./ d; |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
557 sn = norm (s); |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
558 if (sn > 0) |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
559 ## Normalize and rescale. |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
560 s = (s / sn) ./ d; |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
561 ## Get the line minimizer in s direction. |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
562 tn = norm (r*s); |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
563 snm = (sn / tn) / tn; |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
564 if (snm < delta) |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
565 ## Get the dogleg path minimizer. |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
566 bn = norm (b); |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
567 dxn = delta/xn; snmd = snm/delta; |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
568 t = (bn/sn) * (bn/xn) * snmd; |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
569 t -= dxn * snmd^2 - sqrt ((t-dxn)^2 + (1-dxn^2)*(1-snmd^2)); |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
570 alpha = dxn*(1-snmd^2) / t; |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
571 else |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
572 alpha = 0; |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
573 endif |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
574 else |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
575 alpha = delta / xn; |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
576 snm = 0; |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
577 endif |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
578 ## Form the appropriate convex combination. |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
579 x = alpha * x + ((1-alpha) * min (snm, delta)) * s; |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
580 endif |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
581 |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
582 endfunction |
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
583 |
17338
1c89599167a6
maint: End m-files with 1 blank line.
Rik <rik@octave.org>
parents:
17290
diff
changeset
|
584 |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
585 %!function retval = __f (p) |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
586 %! x = p(1); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
587 %! y = p(2); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
588 %! z = p(3); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
589 %! retval = zeros (3, 1); |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
590 %! retval(1) = sin (x) + y^2 + log (z) - 7; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
591 %! retval(2) = 3*x + 2^y -z^3 + 1; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
592 %! retval(3) = x + y + z - 5; |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
593 %!endfunction |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
594 %!test |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
595 %! x_opt = [ 0.599054; |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
596 %! 2.395931; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
597 %! 2.005014 ]; |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
598 %! tol = 1.0e-5; |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
599 %! [x, fval, info] = fsolve (@__f, [ 0.5; 2.0; 2.5 ]); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
600 %! assert (info > 0); |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
601 %! assert (norm (x - x_opt, Inf) < tol); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
602 %! assert (norm (fval) < tol); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
603 |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
604 %!function retval = __f (p) |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
605 %! x = p(1); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
606 %! y = p(2); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
607 %! z = p(3); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
608 %! w = p(4); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
609 %! retval = zeros (4, 1); |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
610 %! retval(1) = 3*x + 4*y + exp (z + w) - 1.007; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
611 %! retval(2) = 6*x - 4*y + exp (3*z + w) - 11; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
612 %! retval(3) = x^4 - 4*y^2 + 6*z - 8*w - 20; |
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
613 %! retval(4) = x^2 + 2*y^3 + z - w - 4; |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
614 %!endfunction |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
615 %!test |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
616 %! x_opt = [ -0.767297326653401, 0.590671081117440, ... |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
617 %! 1.47190018629642, -1.52719341133957 ]; |
25730
bf89eea6652e
Overhaul fsolve for better performance and Matlab compatibility (bug #41146).
Rik <rik@octave.org>
parents:
25606
diff
changeset
|
618 %! tol = 1.0e-4; |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
619 %! [x, fval, info] = fsolve (@__f, [-1, 1, 2, -1]); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
620 %! assert (info > 0); |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
621 %! assert (norm (x - x_opt, Inf) < tol); |
8306
43795cf108d0
initial implementation of fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
diff
changeset
|
622 %! assert (norm (fval) < tol); |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
623 |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
624 %!function retval = __f (p) |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
625 %! x = p(1); |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
626 %! y = p(2); |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
627 %! z = p(3); |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
628 %! retval = zeros (3, 1); |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
629 %! retval(1) = sin (x) + y^2 + log (z) - 7; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
630 %! retval(2) = 3*x + 2^y -z^3 + 1; |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
631 %! retval(3) = x + y + z - 5; |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
632 %! retval(4) = x*x + y - z*log (z) - 1.36; |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
633 %!endfunction |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
634 %!test |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
635 %! x_opt = [ 0.599054; |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
636 %! 2.395931; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
637 %! 2.005014 ]; |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
638 %! tol = 1.0e-5; |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
639 %! [x, fval, info] = fsolve (@__f, [ 0.5; 2.0; 2.5 ]); |
8383
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
640 %! assert (info > 0); |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
641 %! assert (norm (x - x_opt, Inf) < tol); |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
642 %! assert (norm (fval) < tol); |
a762d9daa700
allow overdetermined systems in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8306
diff
changeset
|
643 |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
644 %!function retval = __f (p) |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
645 %! x = p(1); |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
646 %! y = p(2); |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
647 %! z = p(3); |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
648 %! retval = zeros (3, 1); |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
649 %! retval(1) = sin (x) + y^2 + log (z) - 7; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
650 %! retval(2) = 3*x + 2^y -z^3 + 1; |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
651 %! retval(3) = x + y + z - 5; |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
652 %!endfunction |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
653 %!test |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
654 %! x_opt = [ 0.599054; |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
655 %! 2.395931; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14143
diff
changeset
|
656 %! 2.005014 ]; |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
657 %! tol = 1.0e-5; |
8590 | 658 %! opt = optimset ("Updating", "qrp"); |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
659 %! [x, fval, info] = fsolve (@__f, [ 0.5; 2.0; 2.5 ], opt); |
8395
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
660 %! assert (info > 0); |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
661 %! assert (norm (x - x_opt, Inf) < tol); |
88fd356b0d95
optionally allow pivoting in fsolve
Jaroslav Hajek <highegg@gmail.com>
parents:
8383
diff
changeset
|
662 %! assert (norm (fval) < tol); |
8590 | 663 |
664 %!test | |
665 %! b0 = 3; | |
666 %! a0 = 0.2; | |
667 %! x = 0:.5:5; | |
668 %! noise = 1e-5 * sin (100*x); | |
669 %! y = exp (-a0*x) + b0 + noise; | |
670 %! c_opt = [a0, b0]; | |
671 %! tol = 1e-5; | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
672 %! |
21568
3d60ed163b70
maint: Eliminate bad spacing around '='.
Rik <rik@octave.org>
parents:
21546
diff
changeset
|
673 %! [c, fval, info, output] = fsolve (@(c) (exp(-c(1)*x) + c(2) - y), [0, 0]); |
8590 | 674 %! assert (info > 0); |
675 %! assert (norm (c - c_opt, Inf) < tol); | |
676 %! assert (norm (fval) < norm (noise)); | |
677 | |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30875
diff
changeset
|
678 %!function y = cfcn (x) |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
679 %! y(1) = (1+i)*x(1)^2 - (1-i)*x(2) - 2; |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
680 %! y(2) = sqrt (x(1)*x(2)) - (1-2i)*x(3) + (3-4i); |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
681 %! y(3) = x(1) * x(2) - x(3)^2 + (3+2i); |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13027
diff
changeset
|
682 %!endfunction |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
683 |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
684 %!test |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
685 %! x_opt = [-1+i, 1-i, 2+i]; |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
686 %! x = [i, 1, 1+i]; |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
687 %! |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30875
diff
changeset
|
688 %! [x, f, info] = fsolve (@cfcn, x, optimset ("ComplexEqn", "on")); |
8693
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
689 %! tol = 1e-5; |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
690 %! assert (norm (f) < tol); |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
691 %! assert (norm (x - x_opt, Inf) < tol); |
e5ffb52c9c61
improve fsolve and add ComplexEqn option
Jaroslav Hajek <highegg@gmail.com>
parents:
8661
diff
changeset
|
692 |
25499
c29a2107c559
fsolve.m: Add BIST test for bug #53991.
Rik <rik@octave.org>
parents:
25486
diff
changeset
|
693 %!test <*53991> |
c29a2107c559
fsolve.m: Add BIST test for bug #53991.
Rik <rik@octave.org>
parents:
25486
diff
changeset
|
694 %! A = @(lam) [0 1 0 0; 0 0 1 0; 0 0 0 1; 0 0 -lam^2 0]; |
c29a2107c559
fsolve.m: Add BIST test for bug #53991.
Rik <rik@octave.org>
parents:
25486
diff
changeset
|
695 %! C = [1 0 0 0; 0 0 1 0]; |
c29a2107c559
fsolve.m: Add BIST test for bug #53991.
Rik <rik@octave.org>
parents:
25486
diff
changeset
|
696 %! B = @(lam) [C*expm(A(lam)*0); C*expm(A(lam)*1)]; |
c29a2107c559
fsolve.m: Add BIST test for bug #53991.
Rik <rik@octave.org>
parents:
25486
diff
changeset
|
697 %! detB = @(lam) det (B(lam)); |
25606
959657cab846
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
25500
diff
changeset
|
698 %! |
26075
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
699 %! [x, fval, info] = fsolve (detB, 0); |
25499
c29a2107c559
fsolve.m: Add BIST test for bug #53991.
Rik <rik@octave.org>
parents:
25486
diff
changeset
|
700 %! assert (x == 0); |
26075
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
701 %! assert (fval == -1); |
25499
c29a2107c559
fsolve.m: Add BIST test for bug #53991.
Rik <rik@octave.org>
parents:
25486
diff
changeset
|
702 %! assert (info == -2); |
c29a2107c559
fsolve.m: Add BIST test for bug #53991.
Rik <rik@octave.org>
parents:
25486
diff
changeset
|
703 |
25500
b3c35a130f94
fsolve.m: Return info=1 when initial guess (x0) is correct (bug #53991).
Rik <rik@octave.org>
parents:
25499
diff
changeset
|
704 %!test <*53991> |
26075
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
705 %! [x, fval, info] = fsolve (@(x) 5*x, 0); |
25500
b3c35a130f94
fsolve.m: Return info=1 when initial guess (x0) is correct (bug #53991).
Rik <rik@octave.org>
parents:
25499
diff
changeset
|
706 %! assert (x == 0); |
26075
5d28fc32a7c7
doc: Consistent use of "fval" for optimization functions.
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25828
diff
changeset
|
707 %! assert (fval == 0); |
25500
b3c35a130f94
fsolve.m: Return info=1 when initial guess (x0) is correct (bug #53991).
Rik <rik@octave.org>
parents:
25499
diff
changeset
|
708 %! assert (info == 1); |