annotate scripts/optimization/sqp.m @ 27923:bd51beb6205e

update formatting of copyright notices * Use <https://octave.org/copyright/> instead of <https://octave.org/COPYRIGHT.html/>. * For consistency with other comments in the Octave sources, use C++-style comments for copyright blocks in C and C++ files. * Use delimiters above and below copyright blocks that are appropriate for the language used in the file. * Eliminate extra spacing inside copyright blocks. * lex.ll (looks_like_copyright): Also allow newlines and carriage returns before the word "Copyright". * scripts/mk-doc.pl (gethelp): Also skip empty comment lines. * bp-table.cc, type.m: Adjust tests.
author John W. Eaton <jwe@octave.org>
date Wed, 08 Jan 2020 11:59:41 -0500
parents 1891570abac8
children 28de41192f3c 0a5b15007766
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
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 ##
27919
1891570abac8 update Octave Project Developers copyright for the new year
John W. Eaton <jwe@octave.org>
parents: 27918
diff changeset
3 ## Copyright (C) 2005-2020 The Octave Project Developers
27918
b442ec6dda5c use centralized file for copyright info for individual contributors
John W. Eaton <jwe@octave.org>
parents: 26619
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/>.
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
7 ##
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
8 ## This file is part of Octave.
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
9 ##
24534
194eb4bd202b maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents: 23220
diff changeset
10 ## Octave is free software: you can redistribute it and/or modify it
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
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: 23220
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.
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
14 ##
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
15 ## Octave is distributed in the hope that it will be useful, but
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
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.
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
19 ##
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
20 ## You should have received a copy of the GNU General Public License
7016
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 7001
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: 23220
diff changeset
22 ## <https://www.gnu.org/licenses/>.
27923
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
23 ##
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
24 ########################################################################
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
25
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
26 ## -*- texinfo -*-
20852
516bb87ea72e 2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents: 20735
diff changeset
27 ## @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
28 ## @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
29 ## @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
30 ## @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
31 ## @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
32 ## @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
33 ## 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
34 ##
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
35 ## Solve the nonlinear program
6741
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
36 ## @tex
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
37 ## $$
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
38 ## \min_x \phi (x)
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
39 ## $$
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
40 ## @end tex
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
41 ## @ifnottex
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
42 ##
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
43 ## @example
9051
1bf0ce0930be Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents: 8920
diff changeset
44 ## @group
14327
4d917a6a858b doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
45 ## min phi (x)
4d917a6a858b doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
46 ## x
9051
1bf0ce0930be Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents: 8920
diff changeset
47 ## @end group
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
48 ## @end example
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
49 ##
6741
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
50 ## @end ifnottex
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
51 ## subject to
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
52 ## @tex
6741
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
53 ## $$
8167
17352ccd860e describe additional arguments in sqp() documentation string
Ivan Sutoris <ivan.sutoris@gmail.com>
parents: 8047
diff changeset
54 ## g(x) = 0 \qquad h(x) \geq 0 \qquad lb \leq x \leq ub
6741
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
55 ## $$
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
56 ## @end tex
6741
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
57 ## @ifnottex
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
58 ##
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
59 ## @example
9051
1bf0ce0930be Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents: 8920
diff changeset
60 ## @group
14327
4d917a6a858b doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
61 ## g(x) = 0
4d917a6a858b doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
62 ## h(x) >= 0
4d917a6a858b doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
63 ## lb <= x <= ub
9051
1bf0ce0930be Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents: 8920
diff changeset
64 ## @end group
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
65 ## @end example
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10678
diff changeset
66 ##
6741
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
67 ## @end ifnottex
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
68 ## @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
69 ## using a sequential quadratic programming method.
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
70 ##
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
71 ## The first argument is the initial guess for the vector @var{x0}.
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
72 ##
20165
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
73 ## 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
74 ## @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
75 ## return a scalar.
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
76 ##
20165
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
77 ## 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
78 ## 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
79 ## 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
80 ## 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
81 ## 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
82 ## 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
83 ## 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
84 ## Hessian.
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
85 ##
20165
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
86 ## 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
87 ## 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
88 ## @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
89 ##
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
90 ## 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
91 ## 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
92 ## 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
93 ## 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
94 ## @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
95 ## 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
96 ##
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
97 ## The third and fourth arguments may also be 2-element cell arrays of
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
98 ## function handles. The first element should point to the constraint
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
99 ## function and the second should point to a function that computes the
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
100 ## gradient of the constraint function:
6741
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
101 ## @tex
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
102 ## $$
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
103 ## \Bigg( {\partial f(x) \over \partial x_1},
6741
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
104 ## {\partial f(x) \over \partial x_2}, \ldots,
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
105 ## {\partial f(x) \over \partial x_N} \Bigg)^T
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
106 ## $$
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
107 ## @end tex
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
108 ## @ifnottex
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10678
diff changeset
109 ##
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
110 ## @example
9051
1bf0ce0930be Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents: 8920
diff changeset
111 ## @group
14327
4d917a6a858b doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
112 ## [ 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
113 ## transpose ( [ ------ ----- ... ------ ] )
4d917a6a858b doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
114 ## [ dx_1 dx_2 dx_N ]
9051
1bf0ce0930be Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents: 8920
diff changeset
115 ## @end group
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
116 ## @end example
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10678
diff changeset
117 ##
6741
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
118 ## @end ifnottex
20165
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
119 ## 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
120 ## 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
121 ## 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
122 ## 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
123 ## 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
124 ## 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
125 ## 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
126 ##
13207
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
127 ## 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
128 ## 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
129 ##
20165
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
130 ## 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
131 ## 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
132 ##
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
133 ## 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
134 ##
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
135 ## @table @asis
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
136 ## @item 101
19597
db92e7e28e1f strip trailing whitespace from most source files
John W. Eaton <jwe@octave.org>
parents: 19596
diff changeset
137 ## The algorithm terminated normally.
18255
f58a6cd3f909 sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents: 17744
diff changeset
138 ## 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
139 ##
f58a6cd3f909 sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents: 17744
diff changeset
140 ## @item 102
f58a6cd3f909 sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents: 17744
diff changeset
141 ## The BFGS update failed.
f58a6cd3f909 sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents: 17744
diff changeset
142 ##
f58a6cd3f909 sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents: 17744
diff changeset
143 ## @item 103
f58a6cd3f909 sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents: 17744
diff changeset
144 ## 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
145 ##
f58a6cd3f909 sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents: 17744
diff changeset
146 ## @item 104
19597
db92e7e28e1f strip trailing whitespace from most source files
John W. Eaton <jwe@octave.org>
parents: 19596
diff changeset
147 ## 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
148 ## @tex
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
149 ## $\Delta x,$
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
150 ## @end tex
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
151 ## @ifnottex
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
152 ## delta @var{x},
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
153 ## @end ifnottex
13207
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
154 ## 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
155 ## @end table
8167
17352ccd860e describe additional arguments in sqp() documentation string
Ivan Sutoris <ivan.sutoris@gmail.com>
parents: 8047
diff changeset
156 ##
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
157 ## An example of calling @code{sqp}:
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
158 ##
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
159 ## @example
7031
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
160 ## function r = g (x)
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
161 ## r = [ sumsq(x)-10;
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
162 ## x(2)*x(3)-5*x(4)*x(5);
7031
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
163 ## x(1)^3+x(2)^3+1 ];
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
164 ## endfunction
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
165 ##
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
166 ## 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
167 ## obj = exp (prod (x)) - 0.5*(x(1)^3+x(2)^3+1)^2;
7031
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
168 ## endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
169 ##
7031
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
170 ## x0 = [-1.8; 1.7; 1.9; -0.8; -0.8];
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
171 ##
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
172 ## [x, obj, info, iter, nf, lambda] = sqp (x0, @@phi, @@g, [])
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
173 ##
7031
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
174 ## x =
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
175 ##
7031
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
176 ## -1.71714
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
177 ## 1.59571
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
178 ## 1.82725
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
179 ## -0.76364
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
180 ## -0.76364
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
181 ##
7031
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
182 ## obj = 0.053950
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
183 ## info = 101
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
184 ## iter = 8
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
185 ## nf = 10
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
186 ## lambda =
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
187 ##
7031
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
188 ## -0.0401627
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
189 ## 0.0379578
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
190 ## -0.0052227
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
191 ## @end example
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
192 ##
5642
2618a0750ae6 [project @ 2006-03-06 21:26:48 by jwe]
jwe
parents: 5307
diff changeset
193 ## @seealso{qp}
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
194 ## @end deftypefn
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
195
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
196 function [x, obj, info, iter, nf, lambda] = sqp (x0, objf, cef, cif, lb, ub, maxiter, tolerance)
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
197
17336
b81b9d079515 Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents: 17245
diff changeset
198 globals = struct (); # data and handles, needed and changed by subfunctions
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
199
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
200 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
201 print_usage ();
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
202 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
203
19833
9fc020886ae9 maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19697
diff changeset
204 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
205 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
206 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
207 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
208 x0 = x0';
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
209 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
210
26237
34cc6edf2041 sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents: 25054
diff changeset
211 have_grd = 0;
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
212 have_hess = 0;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
213 if (iscell (objf))
11120
a44f979a35ce style fixes for some .m files
John W. Eaton <jwe@octave.org>
parents: 10821
diff changeset
214 switch (numel (objf))
17174
c3c1ebfaa7dc maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents: 16933
diff changeset
215 case 1
c3c1ebfaa7dc maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents: 16933
diff changeset
216 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
217 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
218 case 2
c3c1ebfaa7dc maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents: 16933
diff changeset
219 obj_fun = objf{1};
c3c1ebfaa7dc maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents: 16933
diff changeset
220 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
221 have_grd = 1;
17174
c3c1ebfaa7dc maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents: 16933
diff changeset
222 case 3
c3c1ebfaa7dc maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents: 16933
diff changeset
223 obj_fun = objf{1};
c3c1ebfaa7dc maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents: 16933
diff changeset
224 obj_grd = objf{2};
c3c1ebfaa7dc maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents: 16933
diff changeset
225 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
226 have_grd = 1;
17174
c3c1ebfaa7dc maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents: 16933
diff changeset
227 have_hess = 1;
c3c1ebfaa7dc maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents: 16933
diff changeset
228 otherwise
c3c1ebfaa7dc maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents: 16933
diff changeset
229 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
230 endswitch
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
231 else
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
232 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
233 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
234 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
235
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
236 ce_fun = @empty_cf;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
237 ce_grd = @empty_jac;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
238 if (nargin > 2)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
239 if (iscell (cef))
11120
a44f979a35ce style fixes for some .m files
John W. Eaton <jwe@octave.org>
parents: 10821
diff changeset
240 switch (numel (cef))
17174
c3c1ebfaa7dc maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents: 16933
diff changeset
241 case 1
c3c1ebfaa7dc maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents: 16933
diff changeset
242 ce_fun = cef{1};
19833
9fc020886ae9 maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19697
diff changeset
243 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
244 case 2
c3c1ebfaa7dc maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents: 16933
diff changeset
245 ce_fun = cef{1};
c3c1ebfaa7dc maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents: 16933
diff changeset
246 ce_grd = cef{2};
c3c1ebfaa7dc maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents: 16933
diff changeset
247 otherwise
c3c1ebfaa7dc maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents: 16933
diff changeset
248 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
249 endswitch
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
250 elseif (! isempty (cef))
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
251 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
252 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
253 endif
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
254 endif
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
255
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
256 ci_fun = @empty_cf;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
257 ci_grd = @empty_jac;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
258 if (nargin > 3)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
259 ## 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
260 globals.cif = cif;
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
261 ## 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
262 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
263 if (iscell (cif))
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
264 if (length (cif) > 0)
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
265 globals.cifcn = cif{1};
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
266 endif
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
267 elseif (! isempty (cif))
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
268 globals.cifcn = cif;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
269 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
270
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
271 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
272 ## 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
273 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
274 if (iscell (cif))
16933
e39f00a32dc7 maint: Use parentheses around condition for switch(),while(),if() statements.
Rik <rik@octave.org>
parents: 16064
diff changeset
275 switch (length (cif))
17174
c3c1ebfaa7dc maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents: 16933
diff changeset
276 case 1
c3c1ebfaa7dc maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents: 16933
diff changeset
277 ci_fun = cif{1};
c3c1ebfaa7dc maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents: 16933
diff changeset
278 case 2
c3c1ebfaa7dc maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents: 16933
diff changeset
279 ci_fun = cif{1};
c3c1ebfaa7dc maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents: 16933
diff changeset
280 ci_grd = cif{2};
c3c1ebfaa7dc maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents: 16933
diff changeset
281 otherwise
c3c1ebfaa7dc maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents: 16933
diff changeset
282 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
283 endswitch
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
284 elseif (! isempty (cif))
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
285 ci_fun = cif; # No cell array, only constraint inequality function set
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
286 endif
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
287 else
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
288 ## 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
289 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
290 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
291 if (isvector (lb))
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
292 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
293 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
294 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
295 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
296 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
297 if (isa (x0, "single"))
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
298 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
299 else
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
300 globals.lb = tmp_lb = -realmax;
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10542
diff changeset
301 endif
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
302 else
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
303 error ("sqp: invalid lower bound");
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
304 endif
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
305
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
306 if (isvector (ub))
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
307 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
308 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
309 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
310 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
311 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
312 if (isa (x0, "single"))
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
313 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
314 else
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
315 globals.ub = tmp_ub = realmax;
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10542
diff changeset
316 endif
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
317 else
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
318 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
319 endif
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
320
11347
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
321 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
322 error ("sqp: upper bound smaller than lower bound");
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
323 endif
11347
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
324 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
325 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
326 ci_grd = @(x) cigrad_ub_lb (x, bounds_grad, globals);
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
327 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
328
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
329 endif # if (nargin > 3)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
330
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
331 iter_max = 100;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
332 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
333 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
334 iter_max = maxiter;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
335 else
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
336 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
337 endif
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
338 endif
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
339
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
340 tol = sqrt (eps);
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
341 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
342 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
343 tol = tolerance;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
344 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
345 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
346 endif
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
347 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
348
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
349 ## 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
350 ## 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
351 ## 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
352 ##
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
353 ## obj_fun -- objective function
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
354 ## obj_grad -- objective gradient
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
355 ## 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
356 ## 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
357 ## 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
358 x = x0;
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
359
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
360 obj = feval (obj_fun, x0);
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
361 globals.nfun = 1;
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
362
26237
34cc6edf2041 sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents: 25054
diff changeset
363 if (have_grd)
34cc6edf2041 sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents: 25054
diff changeset
364 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
365 else
34cc6edf2041 sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents: 25054
diff changeset
366 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
367 endif
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
368
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
369 ## 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
370 n = length (x0);
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
371 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
372 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
373 else
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
374 B = eye (n, n);
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
375 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
376
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
377 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
378 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
379
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
380 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
381 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
382
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
383 A = [F; C];
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 ## 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
386 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
387
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
388 qp_iter = 1;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
389 alpha = 1;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
390
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
391 info = 0;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
392 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
393 # 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
394 # report (iter, qp_iter, alpha, __sqp_nfun__, obj);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
395
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
396 while (++iter < iter_max)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
397
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
398 ## 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
399 ## order necessary conditions.
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
400 nr_f = rows (F);
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
401
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
402 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
403 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
404
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
405 con = [ce; ci];
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
406
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
407 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
408 t1 = norm (ce);
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
409 t2 = all (ci >= 0);
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
410 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
411 t4 = norm (lambda .* con);
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
412
18255
f58a6cd3f909 sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents: 17744
diff changeset
413 ## 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
414 ## 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
415 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
416 info = 101;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
417 break;
6382
f74e71ef6612 [project @ 2007-03-05 21:29:40 by jwe]
jwe
parents: 6046
diff changeset
418 endif
f74e71ef6612 [project @ 2007-03-05 21:29:40 by jwe]
jwe
parents: 6046
diff changeset
419
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
420 ## 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
421 g = -ce;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
422 d = -ci;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
423
22224
20e0c0b8820c Allow sequential quadratic programs with infeasible QPs (bug #36015).
Lachlan Andrew <lachlanbis@gmail.com>
parents: 21758
diff changeset
424 old_lambda = lambda;
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
425 [p, obj_qp, INFO, lambda] = qp (x, B, c, F, g, [], [], d, C,
26332
d4211810202a qp and sqp: Non-fixed tolerance for qp (bug #53506).
Maor Shutman <maorus12@gmail.com>
parents: 26237
diff changeset
426 Inf (size (d)), struct ("TolX", tol));
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
427
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
428 info = INFO.info;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
429
18255
f58a6cd3f909 sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents: 17744
diff changeset
430 ## 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
431 ## 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
432
13206
658aa9fd8665 sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents: 11587
diff changeset
433 id = "Octave:SQP-QP-subproblem";
658aa9fd8665 sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents: 11587
diff changeset
434 switch (info)
658aa9fd8665 sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents: 11587
diff changeset
435 case 2
658aa9fd8665 sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents: 11587
diff changeset
436 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
437 case 3
658aa9fd8665 sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents: 11587
diff changeset
438 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
439 INFO.solveiter);
658aa9fd8665 sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents: 11587
diff changeset
440 case 6
658aa9fd8665 sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents: 11587
diff changeset
441 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
442 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
443 endswitch
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
444
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
445 ## 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
446 ## merit function phi.
17245
7babcdb9bc13 Use ... instead of \ for line continuation marker.
Stefan Mahr <dac922@gmx.de>
parents: 17174
diff changeset
447 [x_new, alpha, obj_new, globals] = ...
7babcdb9bc13 Use ... instead of \ for line continuation marker.
Stefan Mahr <dac922@gmx.de>
parents: 17174
diff changeset
448 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
449 obj, c, globals);
26617
98afb8bbd1f6 maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 26376
diff changeset
450
26237
34cc6edf2041 sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents: 25054
diff changeset
451 delx = x_new - x;
34cc6edf2041 sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents: 25054
diff changeset
452
34cc6edf2041 sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents: 25054
diff changeset
453 ## 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
454 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
455 info = 104;
34cc6edf2041 sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents: 25054
diff changeset
456 break;
34cc6edf2041 sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents: 25054
diff changeset
457 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
458
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
459 ## 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
460 if (have_grd)
34cc6edf2041 sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents: 25054
diff changeset
461 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
462 else
34cc6edf2041 sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents: 25054
diff changeset
463 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
464 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
465
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
466 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
467 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
468
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
469 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
470 C_new = feval (ci_grd, x_new);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
471
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
472 A_new = [F_new; C_new];
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
473
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
474 ## Set
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
475 ##
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
476 ## s = alpha * p
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
477 ## y = grad_x L (x_new, lambda) - grad_x L (x, lambda})
6527
2a04f026ef54 [project @ 2007-04-16 16:54:09 by jwe]
jwe
parents: 6382
diff changeset
478
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
479 y = c_new - c;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
480
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
481 if (! isempty (A))
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
482 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
483 y -= t;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
484 endif
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
485
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
486 if (have_hess)
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 B = feval (obj_hess, x);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
489
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
490 else
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
491 ## 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
492 delxt = delx';
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
493
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
494 ## 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
495 ## 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
496 d1 = delxt*B*delx;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
497
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
498 t1 = 0.2 * d1;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
499 t2 = delxt*y;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
500
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
501 if (t2 < t1)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
502 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
503 else
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
504 theta = 1;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
505 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
506
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
507 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
508
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
509 d2 = delxt*r;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
510
18255
f58a6cd3f909 sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents: 17744
diff changeset
511 ## 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
512 ## 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
513 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
514 info = 102;
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10542
diff changeset
515 break;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
516 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
517
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
518 B = B - B*delx*delxt*B/d1 + r*r'/d2;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
519
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
520 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
521
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
522 x = x_new;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
523
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
524 obj = obj_new;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
525
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
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
527
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
528 ce = ce_new;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
529 F = F_new;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
530
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
531 ci = ci_new;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
532 C = C_new;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
533
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
534 A = A_new;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
535
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
536 # 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
537
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
538 endwhile
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
539
18255
f58a6cd3f909 sqp.m: Fix return values in sqp (bug #32008).
Arun Giridhar <arungiridhar@gmail.com>
parents: 17744
diff changeset
540 ## 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
541 if (iter >= iter_max)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
542 info = 103;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
543 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
544
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
545 nf = globals.nfun;
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
546
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
547 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
548
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
549
17245
7babcdb9bc13 Use ... instead of \ for line continuation marker.
Stefan Mahr <dac922@gmx.de>
parents: 17174
diff changeset
550 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
551 x, mu, globals)
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
552
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
553 ce = feval (ce_fun, x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
554 ci = feval (ci_fun, x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
555
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
556 idx = ci < 0;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
557
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
558 con = [ce; ci(idx)];
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
559
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
560 if (isempty (obj))
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
561 obj = feval (obj_fun, x);
20735
418ae0cb752f Replace ++,-- with in-place operators for performance.
Rik <rik@octave.org>
parents: 20165
diff changeset
562 globals.nfun += 1;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
563 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
564
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
565 merit = obj;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
566 t = norm (con, 1) / mu;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
567
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
568 if (! isempty (t))
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
569 merit += t;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
570 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
571
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
572 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
573
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
574
17245
7babcdb9bc13 Use ... instead of \ for line continuation marker.
Stefan Mahr <dac922@gmx.de>
parents: 17174
diff changeset
575 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
576 linesearch_L1 (x, p, obj_fun, obj_grd, ce_fun, ci_fun, lambda, obj, c, globals)
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
577
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
578 ## Choose parameters
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
579 ##
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
580 ## eta in the range (0, 0.5)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
581 ## tau in the range (0, 1)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
582
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
583 eta = 0.25;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
584 tau = 0.5;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
585
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
586 delta_bar = sqrt (eps);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
587
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
588 if (isempty (lambda))
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
589 mu = 1 / delta_bar;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
590 else
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
591 mu = 1 / (norm (lambda, Inf) + delta_bar);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
592 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
593
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
594 alpha = 1;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
595
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
596 ce = feval (ce_fun, x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
597
17245
7babcdb9bc13 Use ... instead of \ for line continuation marker.
Stefan Mahr <dac922@gmx.de>
parents: 17174
diff changeset
598 [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
599 mu, globals);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
600
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
601 D_phi_x_mu = c' * p;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
602 d = feval (ci_fun, x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
603 ## only those elements of d corresponding
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
604 ## to violated constraints should be included.
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
605 idx = d < 0;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
606 t = - norm ([ce; d(idx)], 1) / mu;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
607 if (! isempty (t))
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
608 D_phi_x_mu += t;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
609 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
610
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
611 while (1)
17245
7babcdb9bc13 Use ... instead of \ for line continuation marker.
Stefan Mahr <dac922@gmx.de>
parents: 17174
diff changeset
612 [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
613 x+alpha*p, mu, globals);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
614 p2 = phi_x_mu+eta*alpha*D_phi_x_mu;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
615 if (p1 > p2)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
616 ## Reset alpha = tau_alpha * alpha for some tau_alpha in the
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
617 ## range (0, tau).
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
618 tau_alpha = 0.9 * tau; # ??
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
619 alpha = tau_alpha * alpha;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
620 else
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
621 break;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
622 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
623 endwhile
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
624
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
625 x_new = x + alpha * p;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
626
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
627 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
628
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
629
26237
34cc6edf2041 sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents: 25054
diff changeset
630 function grd = fdgrd (f, x, y0)
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
631
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
632 if (! isempty (f))
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
633 nx = length (x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
634 grd = zeros (nx, 1);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
635 deltax = sqrt (eps);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
636 for i = 1:nx
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
637 t = x(i);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
638 x(i) += deltax;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
639 grd(i) = (feval (f, x) - y0) / deltax;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
640 x(i) = t;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
641 endfor
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
642 else
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
643 grd = zeros (0, 1);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
644 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
645
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
646 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
647
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
648
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
649 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
650
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
651 nx = length (x);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
652 if (! isempty (f))
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
653 y0 = feval (f, x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
654 nf = length (y0);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
655 nx = length (x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
656 jac = zeros (nf, nx);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
657 deltax = sqrt (eps);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
658 for i = 1:nx
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
659 t = x(i);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
660 x(i) += deltax;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
661 jac(:,i) = (feval (f, x) - y0) / deltax;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
662 x(i) = t;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
663 endfor
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
664 else
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
665 jac = zeros (0, nx);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
666 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
667
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
668 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
669
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
670
26237
34cc6edf2041 sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents: 25054
diff changeset
671 function grd = fd_obj_grd (x, obj, obj_fun)
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
672
26237
34cc6edf2041 sqp.m: Eliminate redundant objective function calls (patch #8073).
Ken Marek <gm.kmarek@gmail.com>
parents: 25054
diff changeset
673 grd = fdgrd (obj_fun, x, obj);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
674
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
675 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
676
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
677
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
678 function res = empty_cf (x)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
679
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
680 res = zeros (0, 1);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
681
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
682 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
683
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
684
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
685 function res = empty_jac (x)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
686
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
687 res = zeros (0, length (x));
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
688
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
689 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
690
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
691
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
692 function jac = fd_ce_jac (x, ce_fun)
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
693
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
694 jac = fdjac (ce_fun, x);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
695
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
696 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
697
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
698
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
699 function jac = fd_ci_jac (x, cifcn)
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
700
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
701 ## 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
702 jac = fdjac (cifcn, x);
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
703
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
704 endfunction
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
705
7017
a1dbe9d80eee [project @ 2007-10-12 21:27:11 by jwe]
jwe
parents: 7016
diff changeset
706
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
707 function res = cf_ub_lb (x, lbidx, ubidx, globals)
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
708
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
709 ## 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
710 if (isempty (globals.cifcn))
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
711 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
712 else
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
713 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
714 globals.ub-x(ubidx,1)];
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
715 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
716
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
717 endfunction
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
718
7017
a1dbe9d80eee [project @ 2007-10-12 21:27:11 by jwe]
jwe
parents: 7016
diff changeset
719
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
720 function res = cigrad_ub_lb (x, bgrad, globals)
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
721
19833
9fc020886ae9 maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19697
diff changeset
722 cigradfcn = @(x) fd_ci_jac (x, globals.cifcn);
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
723
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
724 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
725 cigradfcn = globals.cif{2};
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
726 endif
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
727
11347
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
728 if (isempty (cigradfcn))
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
729 res = bgrad;
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
730 else
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
731 res = [feval(cigradfcn,x); bgrad];
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
732 endif
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
733
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
734 endfunction
7361
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
735
19833
9fc020886ae9 maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19697
diff changeset
736 ## 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
737 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
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 if (nargin == 0)
21634
96518f623c91 Backed out changeset dcf8922b724b
Mike Miller <mtmiller@octave.org>
parents: 21633
diff changeset
740 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
741 else
21634
96518f623c91 Backed out changeset dcf8922b724b
Mike Miller <mtmiller@octave.org>
parents: 21633
diff changeset
742 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
743 endif
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
744
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
745 endfunction
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
746
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14327
diff changeset
747
19833
9fc020886ae9 maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19697
diff changeset
748 ################################################################################
9fc020886ae9 maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19697
diff changeset
749 ## 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
750
13305
63463570d9fe Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents: 13268
diff changeset
751 %!function r = __g (x)
7371
fe9a44d753d6 [project @ 2008-01-14 19:34:22 by jwe]
jwe
parents: 7361
diff changeset
752 %! r = [sumsq(x)-10;
fe9a44d753d6 [project @ 2008-01-14 19:34:22 by jwe]
jwe
parents: 7361
diff changeset
753 %! x(2)*x(3)-5*x(4)*x(5);
fe9a44d753d6 [project @ 2008-01-14 19:34:22 by jwe]
jwe
parents: 7361
diff changeset
754 %! 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
755 %!endfunction
7361
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
756 %!
13305
63463570d9fe Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents: 13268
diff changeset
757 %!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
758 %! 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
759 %!endfunction
7361
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
760 %!
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
761 %!test
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
762 %!
7361
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
763 %! x0 = [-1.8; 1.7; 1.9; -0.8; -0.8];
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
764 %!
13305
63463570d9fe Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents: 13268
diff changeset
765 %! [x, obj, info, iter, nf, lambda] = sqp (x0, @__phi, @__g, []);
7361
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
766 %!
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
767 %! x_opt = [-1.717143501952599;
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
768 %! 1.595709610928535;
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
769 %! 1.827245880097156;
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
770 %! -0.763643103133572;
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
771 %! -0.763643068453300];
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
772 %!
7371
fe9a44d753d6 [project @ 2008-01-14 19:34:22 by jwe]
jwe
parents: 7361
diff changeset
773 %! obj_opt = 0.0539498477702739;
7361
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
774 %!
15757
534a2c881f45 Increase test tolerance on sqp (bug #37742)
Mike Miller <mtmiller@ieee.org>
parents: 14868
diff changeset
775 %! 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
776 %! 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
777
19833
9fc020886ae9 maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19697
diff changeset
778 ## Test input validation
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
779 %!error sqp ()
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
780 %!error sqp (1)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
781 %!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
782 %!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
783 %!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
784 %!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
785 %!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
786 %!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
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),[], ones (2,2))
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14327
diff changeset
789 %!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
790 %!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
791 %!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
792 %!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
793 %!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
794 %!error sqp (1, cell (3,1), cell (2,1), cell (2,1),[],[],[],-1)