Mercurial > octave
annotate scripts/optimization/sqp.m @ 26617:98afb8bbd1f6 stable
maint: Strip trailing whitespace from source files.
* NEWS, GLCanvas.cc, dialog.cc, error.cc, ft-text-renderer.cc,
sysdep.cc, utils.cc, pt-jit.cc, oct-rl-hist.c, makeUniqueStrings.m,
randi.m, waitbar.m, clearvars.m, sqp.m, fplot.m, openfig.m,
__print_parse_opts__.m, octave-svgconvert.cc: Strip trailing whitespace.
author | John W. Eaton <jwe@octave.org> |
---|---|
date | Thu, 24 Jan 2019 19:18:51 +0000 |
parents | 00f796120a6d |
children | c9e5a09af7af |
rev | line source |
---|---|
26376
00f796120a6d
maint: Update copyright dates in all source files.
John W. Eaton <jwe@octave.org>
parents:
26237
diff
changeset
|
1 ## Copyright (C) 2005-2019 John W. Eaton |
00f796120a6d
maint: Update copyright dates in all source files.
John W. Eaton <jwe@octave.org>
parents:
26237
diff
changeset
|
2 ## Copyright (C) 2013-2019 Arun Giridhar |
5289 | 3 ## |
4 ## This file is part of Octave. | |
5 ## | |
24534
194eb4bd202b
maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents:
23220
diff
changeset
|
6 ## Octave is free software: you can redistribute it and/or modify it |
5289 | 7 ## 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:
23220
diff
changeset
|
8 ## 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
|
9 ## (at your option) any later version. |
5289 | 10 ## |
11 ## Octave is distributed in the hope that it will be useful, but | |
12 ## WITHOUT ANY WARRANTY; without even the implied warranty of | |
22755
3a2b891d0b33
maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents:
22323
diff
changeset
|
13 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
3a2b891d0b33
maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents:
22323
diff
changeset
|
14 ## GNU General Public License for more details. |
5289 | 15 ## |
16 ## You should have received a copy of the GNU General Public License | |
7016 | 17 ## 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:
23220
diff
changeset
|
18 ## <https://www.gnu.org/licenses/>. |
5289 | 19 |
20 ## -*- texinfo -*- | |
20852
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20735
diff
changeset
|
21 ## @deftypefn {} {[@var{x}, @var{obj}, @var{info}, @var{iter}, @var{nf}, @var{lambda}] =} sqp (@var{x0}, @var{phi}) |
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20735
diff
changeset
|
22 ## @deftypefnx {} {[@dots{}] =} sqp (@var{x0}, @var{phi}, @var{g}) |
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20735
diff
changeset
|
23 ## @deftypefnx {} {[@dots{}] =} sqp (@var{x0}, @var{phi}, @var{g}, @var{h}) |
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20735
diff
changeset
|
24 ## @deftypefnx {} {[@dots{}] =} sqp (@var{x0}, @var{phi}, @var{g}, @var{h}, @var{lb}, @var{ub}) |
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20735
diff
changeset
|
25 ## @deftypefnx {} {[@dots{}] =} sqp (@var{x0}, @var{phi}, @var{g}, @var{h}, @var{lb}, @var{ub}, @var{maxiter}) |
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20735
diff
changeset
|
26 ## @deftypefnx {} {[@dots{}] =} sqp (@var{x0}, @var{phi}, @var{g}, @var{h}, @var{lb}, @var{ub}, @var{maxiter}, @var{tol}) |
20165
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
27 ## Minimize an objective function using sequential quadratic programming (SQP). |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
28 ## |
5289 | 29 ## Solve the nonlinear program |
6741 | 30 ## @tex |
31 ## $$ | |
32 ## \min_x \phi (x) | |
33 ## $$ | |
34 ## @end tex | |
35 ## @ifnottex | |
5289 | 36 ## |
37 ## @example | |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
38 ## @group |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
39 ## min phi (x) |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
40 ## x |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
41 ## @end group |
5289 | 42 ## @end example |
43 ## | |
6741 | 44 ## @end ifnottex |
45 ## subject to | |
5289 | 46 ## @tex |
6741 | 47 ## $$ |
8167
17352ccd860e
describe additional arguments in sqp() documentation string
Ivan Sutoris <ivan.sutoris@gmail.com>
parents:
8047
diff
changeset
|
48 ## g(x) = 0 \qquad h(x) \geq 0 \qquad lb \leq x \leq ub |
6741 | 49 ## $$ |
5289 | 50 ## @end tex |
6741 | 51 ## @ifnottex |
5289 | 52 ## |
53 ## @example | |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
54 ## @group |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
55 ## g(x) = 0 |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
56 ## h(x) >= 0 |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
57 ## lb <= x <= ub |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
58 ## @end group |
5289 | 59 ## @end example |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10678
diff
changeset
|
60 ## |
6741 | 61 ## @end ifnottex |
5289 | 62 ## @noindent |
11348
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
63 ## using a sequential quadratic programming method. |
5289 | 64 ## |
11348
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
65 ## The first argument is the initial guess for the vector @var{x0}. |
5289 | 66 ## |
20165
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
67 ## The second argument is a function handle pointing to the objective function |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
68 ## @var{phi}. The objective function must accept one vector argument and |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
69 ## return a scalar. |
5289 | 70 ## |
20165
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
71 ## The second argument may also be a 2- or 3-element cell array of function |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
72 ## handles. The first element should point to the objective function, the |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
73 ## second should point to a function that computes the gradient of the |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
74 ## objective function, and the third should point to a function that computes |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
75 ## the Hessian of the objective function. If the gradient function is not |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
76 ## supplied, the gradient is computed by finite differences. If the Hessian |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
77 ## function is not supplied, a BFGS update formula is used to approximate the |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
78 ## Hessian. |
5289 | 79 ## |
20165
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
80 ## When supplied, the gradient function @code{@var{phi}@{2@}} must accept one |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
81 ## vector argument and return a vector. When supplied, the Hessian function |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
82 ## @code{@var{phi}@{3@}} must accept one vector argument and return a matrix. |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
83 ## |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
84 ## The third and fourth arguments @var{g} and @var{h} are function handles |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
85 ## pointing to functions that compute the equality constraints and the |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
86 ## inequality constraints, respectively. If the problem does not have |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
87 ## equality (or inequality) constraints, then use an empty matrix ([]) for |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
88 ## @var{g} (or @var{h}). When supplied, these equality and inequality |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
89 ## constraint functions must accept one vector argument and return a vector. |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
90 ## |
5289 | 91 ## The third and fourth arguments may also be 2-element cell arrays of |
92 ## function handles. The first element should point to the constraint | |
93 ## function and the second should point to a function that computes the | |
94 ## gradient of the constraint function: | |
6741 | 95 ## @tex |
96 ## $$ | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
97 ## \Bigg( {\partial f(x) \over \partial x_1}, |
6741 | 98 ## {\partial f(x) \over \partial x_2}, \ldots, |
99 ## {\partial f(x) \over \partial x_N} \Bigg)^T | |
100 ## $$ | |
101 ## @end tex | |
102 ## @ifnottex | |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10678
diff
changeset
|
103 ## |
5289 | 104 ## @example |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
105 ## @group |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
106 ## [ d f(x) d f(x) d f(x) ] |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
107 ## transpose ( [ ------ ----- ... ------ ] ) |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
108 ## [ dx_1 dx_2 dx_N ] |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
109 ## @end group |
5289 | 110 ## @end example |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10678
diff
changeset
|
111 ## |
6741 | 112 ## @end ifnottex |
20165
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
113 ## The fifth and sixth arguments, @var{lb} and @var{ub}, contain lower and |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
114 ## upper bounds on @var{x}. These must be consistent with the equality and |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
115 ## inequality constraints @var{g} and @var{h}. If the arguments are vectors |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
116 ## then @var{x}(i) is bound by @var{lb}(i) and @var{ub}(i). A bound can also |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
117 ## be a scalar in which case all elements of @var{x} will share the same |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
118 ## bound. If only one bound (lb, ub) is specified then the other will |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
119 ## default to (-@var{realmax}, +@var{realmax}). |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
120 ## |
13207
88bd1d1d6657
Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
11587
diff
changeset
|
121 ## The seventh argument @var{maxiter} specifies the maximum number of |
13929
9cae456085c2
Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
13305
diff
changeset
|
122 ## iterations. The default value is 100. |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
123 ## |
20165
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
124 ## The eighth argument @var{tol} specifies the tolerance for the stopping |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
125 ## criteria. The default value is @code{sqrt (eps)}. |
8167
17352ccd860e
describe additional arguments in sqp() documentation string
Ivan Sutoris <ivan.sutoris@gmail.com>
parents:
8047
diff
changeset
|
126 ## |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
127 ## The value returned in @var{info} may be one of the following: |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10678
diff
changeset
|
128 ## |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
129 ## @table @asis |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
130 ## @item 101 |
19597
db92e7e28e1f
strip trailing whitespace from most source files
John W. Eaton <jwe@octave.org>
parents:
19596
diff
changeset
|
131 ## The algorithm terminated normally. |
18255
f58a6cd3f909
sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents:
17744
diff
changeset
|
132 ## All constraints meet the specified tolerance. |
f58a6cd3f909
sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents:
17744
diff
changeset
|
133 ## |
f58a6cd3f909
sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents:
17744
diff
changeset
|
134 ## @item 102 |
f58a6cd3f909
sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents:
17744
diff
changeset
|
135 ## The BFGS update failed. |
f58a6cd3f909
sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents:
17744
diff
changeset
|
136 ## |
f58a6cd3f909
sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents:
17744
diff
changeset
|
137 ## @item 103 |
f58a6cd3f909
sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents:
17744
diff
changeset
|
138 ## The maximum number of iterations was reached. |
f58a6cd3f909
sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents:
17744
diff
changeset
|
139 ## |
f58a6cd3f909
sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents:
17744
diff
changeset
|
140 ## @item 104 |
19597
db92e7e28e1f
strip trailing whitespace from most source files
John W. Eaton <jwe@octave.org>
parents:
19596
diff
changeset
|
141 ## The stepsize has become too small, i.e., |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
142 ## @tex |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
143 ## $\Delta x,$ |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
144 ## @end tex |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
145 ## @ifnottex |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
146 ## delta @var{x}, |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
147 ## @end ifnottex |
13207
88bd1d1d6657
Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
11587
diff
changeset
|
148 ## is less than @code{@var{tol} * norm (x)}. |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
149 ## @end table |
8167
17352ccd860e
describe additional arguments in sqp() documentation string
Ivan Sutoris <ivan.sutoris@gmail.com>
parents:
8047
diff
changeset
|
150 ## |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
151 ## An example of calling @code{sqp}: |
5289 | 152 ## |
153 ## @example | |
7031 | 154 ## function r = g (x) |
155 ## r = [ sumsq(x)-10; | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
156 ## x(2)*x(3)-5*x(4)*x(5); |
7031 | 157 ## x(1)^3+x(2)^3+1 ]; |
158 ## endfunction | |
159 ## | |
160 ## function obj = phi (x) | |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
161 ## obj = exp (prod (x)) - 0.5*(x(1)^3+x(2)^3+1)^2; |
7031 | 162 ## endfunction |
5289 | 163 ## |
7031 | 164 ## x0 = [-1.8; 1.7; 1.9; -0.8; -0.8]; |
165 ## | |
166 ## [x, obj, info, iter, nf, lambda] = sqp (x0, @@phi, @@g, []) | |
5289 | 167 ## |
7031 | 168 ## x = |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
169 ## |
7031 | 170 ## -1.71714 |
171 ## 1.59571 | |
172 ## 1.82725 | |
173 ## -0.76364 | |
174 ## -0.76364 | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
175 ## |
7031 | 176 ## obj = 0.053950 |
177 ## info = 101 | |
178 ## iter = 8 | |
179 ## nf = 10 | |
180 ## lambda = | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
181 ## |
7031 | 182 ## -0.0401627 |
183 ## 0.0379578 | |
184 ## -0.0052227 | |
5289 | 185 ## @end example |
186 ## | |
5642 | 187 ## @seealso{qp} |
5289 | 188 ## @end deftypefn |
189 | |
11348
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
190 function [x, obj, info, iter, nf, lambda] = sqp (x0, objf, cef, cif, lb, ub, maxiter, tolerance) |
5289 | 191 |
17336
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
17245
diff
changeset
|
192 globals = struct (); # data and handles, needed and changed by subfunctions |
5289 | 193 |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
194 if (nargin < 2 || nargin > 8 || nargin == 5) |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
195 print_usage (); |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
196 endif |
5289 | 197 |
19833
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
198 if (! isvector (x0)) |
11348
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
199 error ("sqp: X0 must be a vector"); |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
200 endif |
11348
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
201 if (rows (x0) == 1) |
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
202 x0 = x0'; |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
203 endif |
5289 | 204 |
26237
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
205 have_grd = 0; |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
206 have_hess = 0; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
207 if (iscell (objf)) |
11120
a44f979a35ce
style fixes for some .m files
John W. Eaton <jwe@octave.org>
parents:
10821
diff
changeset
|
208 switch (numel (objf)) |
17174
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
209 case 1 |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
210 obj_fun = objf{1}; |
26237
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
211 obj_grd = @(x, obj) fd_obj_grd (x, obj, obj_fun); |
17174
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
212 case 2 |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
213 obj_fun = objf{1}; |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
214 obj_grd = objf{2}; |
26237
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
215 have_grd = 1; |
17174
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
216 case 3 |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
217 obj_fun = objf{1}; |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
218 obj_grd = objf{2}; |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
219 obj_hess = objf{3}; |
26237
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
220 have_grd = 1; |
17174
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
221 have_hess = 1; |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
222 otherwise |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
223 error ("sqp: invalid objective function specification"); |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
224 endswitch |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
225 else |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
226 obj_fun = objf; # No cell array, only obj_fun set |
26237
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
227 obj_grd = @(x, obj) fd_obj_grd (x, obj, obj_fun); |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
228 endif |
5289 | 229 |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
230 ce_fun = @empty_cf; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
231 ce_grd = @empty_jac; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
232 if (nargin > 2) |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
233 if (iscell (cef)) |
11120
a44f979a35ce
style fixes for some .m files
John W. Eaton <jwe@octave.org>
parents:
10821
diff
changeset
|
234 switch (numel (cef)) |
17174
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
235 case 1 |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
236 ce_fun = cef{1}; |
19833
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
237 ce_grd = @(x) fd_ce_jac (x, ce_fun); |
17174
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
238 case 2 |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
239 ce_fun = cef{1}; |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
240 ce_grd = cef{2}; |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
241 otherwise |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
242 error ("sqp: invalid equality constraint function specification"); |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
243 endswitch |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
244 elseif (! isempty (cef)) |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
245 ce_fun = cef; # No cell array, only constraint equality function set |
19833
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
246 ce_grd = @(x) fd_ce_jac (x, ce_fun); |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
247 endif |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
248 endif |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
249 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
250 ci_fun = @empty_cf; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
251 ci_grd = @empty_jac; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
252 if (nargin > 3) |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
253 ## constraint function given by user with possible gradient |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
254 globals.cif = cif; |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
255 ## constraint function given by user without gradient |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
256 globals.cifcn = @empty_cf; |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
257 if (iscell (cif)) |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
258 if (length (cif) > 0) |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
259 globals.cifcn = cif{1}; |
5289 | 260 endif |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
261 elseif (! isempty (cif)) |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
262 globals.cifcn = cif; |
5289 | 263 endif |
264 | |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
265 if (nargin < 5 || (nargin > 5 && isempty (lb) && isempty (ub))) |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
266 ## constraint inequality function only without any bounds |
19833
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
267 ci_grd = @(x) fd_ci_jac (x, globals.cifcn); |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
268 if (iscell (cif)) |
16933
e39f00a32dc7
maint: Use parentheses around condition for switch(),while(),if() statements.
Rik <rik@octave.org>
parents:
16064
diff
changeset
|
269 switch (length (cif)) |
17174
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
270 case 1 |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
271 ci_fun = cif{1}; |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
272 case 2 |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
273 ci_fun = cif{1}; |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
274 ci_grd = cif{2}; |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
275 otherwise |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
276 error ("sqp: invalid inequality constraint function specification"); |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
277 endswitch |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
278 elseif (! isempty (cif)) |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
279 ci_fun = cif; # No cell array, only constraint inequality function set |
5289 | 280 endif |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
281 else |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
282 ## constraint inequality function with bounds present |
11348
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
283 lb_idx = ub_idx = true (size (x0)); |
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
284 ub_grad = - (lb_grad = eye (rows (x0))); |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
285 if (isvector (lb)) |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
286 globals.lb = tmp_lb = lb(:); |
11347
2726132f77f6
sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents:
11290
diff
changeset
|
287 lb_idx(:) = tmp_idx = (lb != -Inf); |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
288 globals.lb = globals.lb(tmp_idx, 1); |
11347
2726132f77f6
sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents:
11290
diff
changeset
|
289 lb_grad = lb_grad(lb_idx, :); |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
290 elseif (isempty (lb)) |
11348
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
291 if (isa (x0, "single")) |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
292 globals.lb = tmp_lb = -realmax ("single"); |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
293 else |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
294 globals.lb = tmp_lb = -realmax; |
10549 | 295 endif |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
296 else |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
297 error ("sqp: invalid lower bound"); |
6768 | 298 endif |
299 | |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
300 if (isvector (ub)) |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
301 globals.ub = tmp_ub = ub(:); |
11347
2726132f77f6
sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents:
11290
diff
changeset
|
302 ub_idx(:) = tmp_idx = (ub != Inf); |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
303 globals.ub = globals.ub(tmp_idx, 1); |
11347
2726132f77f6
sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents:
11290
diff
changeset
|
304 ub_grad = ub_grad(ub_idx, :); |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
305 elseif (isempty (ub)) |
11348
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
306 if (isa (x0, "single")) |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
307 globals.ub = tmp_ub = realmax ("single"); |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
308 else |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
309 globals.ub = tmp_ub = realmax; |
10549 | 310 endif |
6768 | 311 else |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
312 error ("sqp: invalid upper bound"); |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
313 endif |
6768 | 314 |
11347
2726132f77f6
sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents:
11290
diff
changeset
|
315 if (any (tmp_lb > tmp_ub)) |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
316 error ("sqp: upper bound smaller than lower bound"); |
6768 | 317 endif |
11347
2726132f77f6
sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents:
11290
diff
changeset
|
318 bounds_grad = [lb_grad; ub_grad]; |
19833
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
319 ci_fun = @(x) cf_ub_lb (x, lb_idx, ub_idx, globals); |
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
320 ci_grd = @(x) cigrad_ub_lb (x, bounds_grad, globals); |
6768 | 321 endif |
5289 | 322 |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
323 endif # if (nargin > 3) |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
324 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
325 iter_max = 100; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
326 if (nargin > 6 && ! isempty (maxiter)) |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
327 if (isscalar (maxiter) && maxiter > 0 && fix (maxiter) == maxiter) |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
328 iter_max = maxiter; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
329 else |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
330 error ("sqp: invalid number of maximum iterations"); |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
331 endif |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
332 endif |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
333 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
334 tol = sqrt (eps); |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
335 if (nargin > 7 && ! isempty (tolerance)) |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
336 if (isscalar (tolerance) && tolerance > 0) |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
337 tol = tolerance; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
338 else |
11472
1740012184f9
Use uppercase for variable names in error() strings to match Info documentation. Only m-files done.
Rik <octave@nomad.inbox5.com>
parents:
11371
diff
changeset
|
339 error ("sqp: invalid value for TOLERANCE"); |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
340 endif |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
341 endif |
5289 | 342 |
11348
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
343 ## Initialize variables for search loop |
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
344 ## Seed x with initial guess and evaluate objective function, constraints, |
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
345 ## and gradients at initial value x0. |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
346 ## |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
347 ## obj_fun -- objective function |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
348 ## obj_grad -- objective gradient |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
349 ## ce_fun -- equality constraint functions |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
350 ## ci_fun -- inequality constraint functions |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
351 ## A == [grad_{x_1} cx_fun, grad_{x_2} cx_fun, ..., grad_{x_n} cx_fun]^T |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
352 x = x0; |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
353 |
11348
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
354 obj = feval (obj_fun, x0); |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
355 globals.nfun = 1; |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
356 |
26237
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
357 if (have_grd) |
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
358 c = feval (obj_grd, x0); |
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
359 else |
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
360 c = feval (obj_grd, x0, obj); |
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
361 endif |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
362 |
11348
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
363 ## Choose an initial NxN symmetric positive definite Hessian approximation B. |
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
364 n = length (x0); |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
365 if (have_hess) |
11348
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
366 B = feval (obj_hess, x0); |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
367 else |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
368 B = eye (n, n); |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
369 endif |
5289 | 370 |
11348
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
371 ce = feval (ce_fun, x0); |
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
372 F = feval (ce_grd, x0); |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
373 |
11348
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
374 ci = feval (ci_fun, x0); |
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
375 C = feval (ci_grd, x0); |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
376 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
377 A = [F; C]; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
378 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
379 ## Choose an initial lambda (x is provided by the caller). |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
380 lambda = 100 * ones (rows (A), 1); |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
381 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
382 qp_iter = 1; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
383 alpha = 1; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
384 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
385 info = 0; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
386 iter = 0; |
11348
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
387 # report (); # Called with no arguments to initialize reporting |
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
388 # report (iter, qp_iter, alpha, __sqp_nfun__, obj); |
5289 | 389 |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
390 while (++iter < iter_max) |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
391 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
392 ## Check convergence. This is just a simple check on the first |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
393 ## order necessary conditions. |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
394 nr_f = rows (F); |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
395 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
396 lambda_e = lambda((1:nr_f)'); |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
397 lambda_i = lambda((nr_f+1:end)'); |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
398 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
399 con = [ce; ci]; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
400 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
401 t0 = norm (c - A' * lambda); |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
402 t1 = norm (ce); |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
403 t2 = all (ci >= 0); |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
404 t3 = all (lambda_i >= 0); |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
405 t4 = norm (lambda .* con); |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
406 |
18255
f58a6cd3f909
sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents:
17744
diff
changeset
|
407 ## Normal convergence. All constraints are satisfied |
f58a6cd3f909
sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents:
17744
diff
changeset
|
408 ## and objective has converged. |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
409 if (t2 && t3 && max ([t0; t1; t4]) < tol) |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
410 info = 101; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
411 break; |
6382 | 412 endif |
413 | |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
414 ## Compute search direction p by solving QP. |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
415 g = -ce; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
416 d = -ci; |
5289 | 417 |
22224
20e0c0b8820c
Allow sequential quadratic programs with infeasible QPs (bug #36015).
Lachlan Andrew <lachlanbis@gmail.com>
parents:
21758
diff
changeset
|
418 old_lambda = lambda; |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
419 [p, obj_qp, INFO, lambda] = qp (x, B, c, F, g, [], [], d, C, |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
420 Inf (size (d))); |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
421 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
422 info = INFO.info; |
5289 | 423 |
18255
f58a6cd3f909
sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents:
17744
diff
changeset
|
424 ## FIXME: check QP solution and attempt to recover if it has failed. |
f58a6cd3f909
sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents:
17744
diff
changeset
|
425 ## For now, just warn about possible problems. |
19593
446c46af4b42
strip trailing whitespace from most source files
John W. Eaton <jwe@octave.org>
parents:
17744
diff
changeset
|
426 |
13206
658aa9fd8665
sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents:
11587
diff
changeset
|
427 id = "Octave:SQP-QP-subproblem"; |
658aa9fd8665
sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents:
11587
diff
changeset
|
428 switch (info) |
658aa9fd8665
sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents:
11587
diff
changeset
|
429 case 2 |
658aa9fd8665
sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents:
11587
diff
changeset
|
430 warning (id, "sqp: QP subproblem is non-convex and unbounded"); |
658aa9fd8665
sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents:
11587
diff
changeset
|
431 case 3 |
658aa9fd8665
sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents:
11587
diff
changeset
|
432 warning (id, "sqp: QP subproblem failed to converge in %d iterations", |
658aa9fd8665
sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents:
11587
diff
changeset
|
433 INFO.solveiter); |
658aa9fd8665
sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents:
11587
diff
changeset
|
434 case 6 |
658aa9fd8665
sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents:
11587
diff
changeset
|
435 warning (id, "sqp: QP subproblem is infeasible"); |
22224
20e0c0b8820c
Allow sequential quadratic programs with infeasible QPs (bug #36015).
Lachlan Andrew <lachlanbis@gmail.com>
parents:
21758
diff
changeset
|
436 lambda = old_lambda; # The return value was size 0x0 in this case. |
13206
658aa9fd8665
sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents:
11587
diff
changeset
|
437 endswitch |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
438 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
439 ## Choose mu such that p is a descent direction for the chosen |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
440 ## merit function phi. |
17245
7babcdb9bc13
Use ... instead of \ for line continuation marker.
Stefan Mahr <dac922@gmx.de>
parents:
17174
diff
changeset
|
441 [x_new, alpha, obj_new, globals] = ... |
7babcdb9bc13
Use ... instead of \ for line continuation marker.
Stefan Mahr <dac922@gmx.de>
parents:
17174
diff
changeset
|
442 linesearch_L1 (x, p, obj_fun, obj_grd, ce_fun, ci_fun, lambda, ... |
26237
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
443 obj, c, globals); |
26617
98afb8bbd1f6
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
26376
diff
changeset
|
444 |
26237
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
445 delx = x_new - x; |
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
446 |
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
447 ## Check if step size has become too small (indicates lack of progress). |
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
448 if (norm (delx) < tol * norm (x)) |
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
449 info = 104; |
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
450 break; |
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
451 endif |
5289 | 452 |
11348
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
453 ## Evaluate objective function, constraints, and gradients at x_new. |
26237
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
454 if (have_grd) |
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
455 c_new = feval (obj_grd, x_new); |
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
456 else |
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
457 c_new = feval (obj_grd, x_new, obj_new); |
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
458 endif |
5289 | 459 |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
460 ce_new = feval (ce_fun, x_new); |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
461 F_new = feval (ce_grd, x_new); |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
462 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
463 ci_new = feval (ci_fun, x_new); |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
464 C_new = feval (ci_grd, x_new); |
5289 | 465 |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
466 A_new = [F_new; C_new]; |
5289 | 467 |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
468 ## Set |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
469 ## |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
470 ## s = alpha * p |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
471 ## y = grad_x L (x_new, lambda) - grad_x L (x, lambda}) |
6527 | 472 |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
473 y = c_new - c; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
474 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
475 if (! isempty (A)) |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
476 t = ((A_new - A)'*lambda); |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
477 y -= t; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
478 endif |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
479 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
480 if (have_hess) |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
481 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
482 B = feval (obj_hess, x); |
5289 | 483 |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
484 else |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
485 ## Update B using a quasi-Newton formula. |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
486 delxt = delx'; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
487 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
488 ## Damped BFGS. Or maybe we would actually want to use the Hessian |
11348
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
489 ## of the Lagrangian, computed directly? |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
490 d1 = delxt*B*delx; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
491 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
492 t1 = 0.2 * d1; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
493 t2 = delxt*y; |
5289 | 494 |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
495 if (t2 < t1) |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
496 theta = 0.8*d1/(d1 - t2); |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
497 else |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
498 theta = 1; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
499 endif |
5289 | 500 |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
501 r = theta*y + (1-theta)*B*delx; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
502 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
503 d2 = delxt*r; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
504 |
18255
f58a6cd3f909
sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents:
17744
diff
changeset
|
505 ## Check if the next BFGS update will work properly. |
f58a6cd3f909
sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents:
17744
diff
changeset
|
506 ## If d1 or d2 vanish, the BFGS update will fail. |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
507 if (d1 == 0 || d2 == 0) |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
508 info = 102; |
10549 | 509 break; |
5289 | 510 endif |
511 | |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
512 B = B - B*delx*delxt*B/d1 + r*r'/d2; |
5289 | 513 |
514 endif | |
515 | |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
516 x = x_new; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
517 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
518 obj = obj_new; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
519 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
520 c = c_new; |
5289 | 521 |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
522 ce = ce_new; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
523 F = F_new; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
524 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
525 ci = ci_new; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
526 C = C_new; |
5289 | 527 |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
528 A = A_new; |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
529 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
530 # report (iter, qp_iter, alpha, __sqp_nfun__, obj); |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
531 |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
532 endwhile |
5289 | 533 |
18255
f58a6cd3f909
sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents:
17744
diff
changeset
|
534 ## Check if we've spent too many iterations without converging. |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
535 if (iter >= iter_max) |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
536 info = 103; |
5289 | 537 endif |
538 | |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
539 nf = globals.nfun; |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
540 |
7399 | 541 endfunction |
5289 | 542 |
543 | |
17245
7babcdb9bc13
Use ... instead of \ for line continuation marker.
Stefan Mahr <dac922@gmx.de>
parents:
17174
diff
changeset
|
544 function [merit, obj, globals] = phi_L1 (obj, obj_fun, ce_fun, ci_fun, ... |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
545 x, mu, globals) |
5289 | 546 |
547 ce = feval (ce_fun, x); | |
548 ci = feval (ci_fun, x); | |
549 | |
550 idx = ci < 0; | |
551 | |
552 con = [ce; ci(idx)]; | |
553 | |
554 if (isempty (obj)) | |
555 obj = feval (obj_fun, x); | |
20735
418ae0cb752f
Replace ++,-- with in-place operators for performance.
Rik <rik@octave.org>
parents:
20165
diff
changeset
|
556 globals.nfun += 1; |
5289 | 557 endif |
558 | |
559 merit = obj; | |
560 t = norm (con, 1) / mu; | |
561 | |
562 if (! isempty (t)) | |
563 merit += t; | |
564 endif | |
565 | |
7399 | 566 endfunction |
5289 | 567 |
568 | |
17245
7babcdb9bc13
Use ... instead of \ for line continuation marker.
Stefan Mahr <dac922@gmx.de>
parents:
17174
diff
changeset
|
569 function [x_new, alpha, obj, globals] = ... |
26237
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
570 linesearch_L1 (x, p, obj_fun, obj_grd, ce_fun, ci_fun, lambda, obj, c, globals) |
5289 | 571 |
572 ## Choose parameters | |
573 ## | |
574 ## eta in the range (0, 0.5) | |
575 ## tau in the range (0, 1) | |
576 | |
577 eta = 0.25; | |
578 tau = 0.5; | |
579 | |
580 delta_bar = sqrt (eps); | |
581 | |
582 if (isempty (lambda)) | |
583 mu = 1 / delta_bar; | |
584 else | |
585 mu = 1 / (norm (lambda, Inf) + delta_bar); | |
586 endif | |
587 | |
588 alpha = 1; | |
589 | |
590 ce = feval (ce_fun, x); | |
591 | |
17245
7babcdb9bc13
Use ... instead of \ for line continuation marker.
Stefan Mahr <dac922@gmx.de>
parents:
17174
diff
changeset
|
592 [phi_x_mu, obj, globals] = phi_L1 (obj, obj_fun, ce_fun, ci_fun, x, ... |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
593 mu, globals); |
5289 | 594 |
595 D_phi_x_mu = c' * p; | |
596 d = feval (ci_fun, x); | |
597 ## only those elements of d corresponding | |
598 ## to violated constraints should be included. | |
599 idx = d < 0; | |
600 t = - norm ([ce; d(idx)], 1) / mu; | |
601 if (! isempty (t)) | |
602 D_phi_x_mu += t; | |
603 endif | |
604 | |
605 while (1) | |
17245
7babcdb9bc13
Use ... instead of \ for line continuation marker.
Stefan Mahr <dac922@gmx.de>
parents:
17174
diff
changeset
|
606 [p1, obj, globals] = phi_L1 ([], obj_fun, ce_fun, ci_fun, ... |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
607 x+alpha*p, mu, globals); |
5289 | 608 p2 = phi_x_mu+eta*alpha*D_phi_x_mu; |
609 if (p1 > p2) | |
610 ## Reset alpha = tau_alpha * alpha for some tau_alpha in the | |
611 ## range (0, tau). | |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
612 tau_alpha = 0.9 * tau; # ?? |
5289 | 613 alpha = tau_alpha * alpha; |
614 else | |
615 break; | |
616 endif | |
617 endwhile | |
618 | |
619 x_new = x + alpha * p; | |
620 | |
7399 | 621 endfunction |
5289 | 622 |
623 | |
26237
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
624 function grd = fdgrd (f, x, y0) |
5289 | 625 |
626 if (! isempty (f)) | |
627 nx = length (x); | |
628 grd = zeros (nx, 1); | |
629 deltax = sqrt (eps); | |
630 for i = 1:nx | |
631 t = x(i); | |
632 x(i) += deltax; | |
633 grd(i) = (feval (f, x) - y0) / deltax; | |
634 x(i) = t; | |
635 endfor | |
636 else | |
637 grd = zeros (0, 1); | |
638 endif | |
639 | |
7399 | 640 endfunction |
5289 | 641 |
642 | |
643 function jac = fdjac (f, x) | |
11348
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
644 |
6768 | 645 nx = length (x); |
5289 | 646 if (! isempty (f)) |
647 y0 = feval (f, x); | |
648 nf = length (y0); | |
649 nx = length (x); | |
650 jac = zeros (nf, nx); | |
651 deltax = sqrt (eps); | |
652 for i = 1:nx | |
653 t = x(i); | |
654 x(i) += deltax; | |
655 jac(:,i) = (feval (f, x) - y0) / deltax; | |
656 x(i) = t; | |
657 endfor | |
658 else | |
659 jac = zeros (0, nx); | |
660 endif | |
661 | |
7399 | 662 endfunction |
5289 | 663 |
664 | |
26237
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
665 function grd = fd_obj_grd (x, obj, obj_fun) |
5289 | 666 |
26237
34cc6edf2041
sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents:
25054
diff
changeset
|
667 grd = fdgrd (obj_fun, x, obj); |
5289 | 668 |
7399 | 669 endfunction |
5289 | 670 |
671 | |
672 function res = empty_cf (x) | |
673 | |
674 res = zeros (0, 1); | |
675 | |
7399 | 676 endfunction |
5289 | 677 |
678 | |
679 function res = empty_jac (x) | |
680 | |
681 res = zeros (0, length (x)); | |
682 | |
7399 | 683 endfunction |
5289 | 684 |
685 | |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
686 function jac = fd_ce_jac (x, ce_fun) |
5289 | 687 |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
688 jac = fdjac (ce_fun, x); |
5289 | 689 |
7399 | 690 endfunction |
5289 | 691 |
692 | |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
693 function jac = fd_ci_jac (x, cifcn) |
5289 | 694 |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
695 ## cifcn = constraint function without gradients and lb or ub |
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
696 jac = fdjac (cifcn, x); |
6768 | 697 |
7399 | 698 endfunction |
6768 | 699 |
7017 | 700 |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
701 function res = cf_ub_lb (x, lbidx, ubidx, globals) |
5289 | 702 |
6768 | 703 ## combine constraint function with ub and lb |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
704 if (isempty (globals.cifcn)) |
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
705 res = [x(lbidx,1)-globals.lb; globals.ub-x(ubidx,1)]; |
11347
2726132f77f6
sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents:
11290
diff
changeset
|
706 else |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
707 res = [feval(globals.cifcn,x); x(lbidx,1)-globals.lb; |
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
708 globals.ub-x(ubidx,1)]; |
6768 | 709 endif |
5289 | 710 |
7399 | 711 endfunction |
6768 | 712 |
7017 | 713 |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
714 function res = cigrad_ub_lb (x, bgrad, globals) |
6768 | 715 |
19833
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
716 cigradfcn = @(x) fd_ci_jac (x, globals.cifcn); |
6768 | 717 |
14595
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
718 if (iscell (globals.cif) && length (globals.cif) > 1) |
11a9d448fdc3
Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents:
14552
diff
changeset
|
719 cigradfcn = globals.cif{2}; |
6768 | 720 endif |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
721 |
11347
2726132f77f6
sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents:
11290
diff
changeset
|
722 if (isempty (cigradfcn)) |
2726132f77f6
sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents:
11290
diff
changeset
|
723 res = bgrad; |
2726132f77f6
sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents:
11290
diff
changeset
|
724 else |
2726132f77f6
sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents:
11290
diff
changeset
|
725 res = [feval(cigradfcn,x); bgrad]; |
6768 | 726 endif |
727 | |
7399 | 728 endfunction |
7361 | 729 |
19833
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
730 ## Utility function used to debug sqp |
11348
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
731 function report (iter, qp_iter, alpha, nfun, obj) |
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
732 |
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
733 if (nargin == 0) |
21634
96518f623c91
Backed out changeset dcf8922b724b
Mike Miller <mtmiller@octave.org>
parents:
21633
diff
changeset
|
734 printf (" Itn ItQP Step Nfun Objective\n"); |
11348
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
735 else |
21634
96518f623c91
Backed out changeset dcf8922b724b
Mike Miller <mtmiller@octave.org>
parents:
21633
diff
changeset
|
736 printf ("%5d %4d %8.1g %5d %13.6e\n", iter, qp_iter, alpha, nfun, obj); |
11348
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
737 endif |
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
738 |
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
739 endfunction |
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
740 |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
741 |
19833
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
742 ################################################################################ |
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
743 ## Test Code |
11348
2ae0ca4ee36b
sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents:
11347
diff
changeset
|
744 |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13268
diff
changeset
|
745 %!function r = __g (x) |
7371 | 746 %! r = [sumsq(x)-10; |
747 %! x(2)*x(3)-5*x(4)*x(5); | |
748 %! x(1)^3+x(2)^3+1 ]; | |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13268
diff
changeset
|
749 %!endfunction |
7361 | 750 %! |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13268
diff
changeset
|
751 %!function obj = __phi (x) |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
752 %! obj = exp (prod (x)) - 0.5*(x(1)^3 + x(2)^3 + 1)^2; |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13268
diff
changeset
|
753 %!endfunction |
7361 | 754 %! |
755 %!test | |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
756 %! |
7361 | 757 %! x0 = [-1.8; 1.7; 1.9; -0.8; -0.8]; |
758 %! | |
13305
63463570d9fe
Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents:
13268
diff
changeset
|
759 %! [x, obj, info, iter, nf, lambda] = sqp (x0, @__phi, @__g, []); |
7361 | 760 %! |
761 %! x_opt = [-1.717143501952599; | |
762 %! 1.595709610928535; | |
763 %! 1.827245880097156; | |
764 %! -0.763643103133572; | |
765 %! -0.763643068453300]; | |
766 %! | |
7371 | 767 %! obj_opt = 0.0539498477702739; |
7361 | 768 %! |
15757
534a2c881f45
Increase test tolerance on sqp (bug #37742)
Mike Miller <mtmiller@ieee.org>
parents:
14868
diff
changeset
|
769 %! assert (x, x_opt, 8*sqrt (eps)); |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
770 %! assert (obj, obj_opt, sqrt (eps)); |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
771 |
19833
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
772 ## Test input validation |
10678
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
773 %!error sqp () |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
774 %!error sqp (1) |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
775 %!error sqp (1,2,3,4,5,6,7,8,9) |
35338deff753
Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
776 %!error sqp (1,2,3,4,5) |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
777 %!error sqp (ones (2,2)) |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
778 %!error sqp (1, cell (4,1)) |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
779 %!error sqp (1, cell (3,1), cell (3,1)) |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
780 %!error sqp (1, cell (3,1), cell (2,1), cell (3,1)) |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
781 %!error sqp (1, cell (3,1), cell (2,1), cell (2,1), ones (2,2),[]) |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
782 %!error sqp (1, cell (3,1), cell (2,1), cell (2,1),[], ones (2,2)) |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
783 %!error sqp (1, cell (3,1), cell (2,1), cell (2,1),1,-1) |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
784 %!error sqp (1, cell (3,1), cell (2,1), cell (2,1),[],[], ones (2,2)) |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
785 %!error sqp (1, cell (3,1), cell (2,1), cell (2,1),[],[],-1) |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
786 %!error sqp (1, cell (3,1), cell (2,1), cell (2,1),[],[],1.5) |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
787 %!error sqp (1, cell (3,1), cell (2,1), cell (2,1),[],[],[], ones (2,2)) |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
788 %!error sqp (1, cell (3,1), cell (2,1), cell (2,1),[],[],[],-1) |