annotate scripts/optimization/sqp.m @ 16064:9bfbe785e726

maint: periodic merge of stable to default
author John W. Eaton <jwe@octave.org>
date Mon, 11 Feb 2013 16:26:08 -0500
parents 534a2c881f45 69c0728def22
children e39f00a32dc7
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
14138
72c96de7a403 maint: update copyright notices for 2012
John W. Eaton <jwe@octave.org>
parents: 13931
diff changeset
1 ## Copyright (C) 2005-2012 John W. Eaton
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
2 ##
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
3 ## This file is part of Octave.
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
4 ##
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
5 ## Octave is free software; you can redistribute it and/or modify it
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
6 ## under the terms of the GNU General Public License as published by
7016
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 7001
diff changeset
7 ## the Free Software Foundation; either version 3 of the License, or (at
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 7001
diff changeset
8 ## your option) any later version.
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
9 ##
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
10 ## 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
11 ## WITHOUT ANY WARRANTY; without even the implied warranty of
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
12 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
13 ## General Public License for more details.
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 ## 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
16 ## along with Octave; see the file COPYING. If not, see
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 7001
diff changeset
17 ## <http://www.gnu.org/licenses/>.
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
18
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
19 ## -*- texinfo -*-
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
20 ## @deftypefn {Function File} {[@var{x}, @var{obj}, @var{info}, @var{iter}, @var{nf}, @var{lambda}] =} sqp (@var{x0}, @var{phi})
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
21 ## @deftypefnx {Function File} {[@dots{}] =} sqp (@var{x0}, @var{phi}, @var{g})
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
22 ## @deftypefnx {Function File} {[@dots{}] =} sqp (@var{x0}, @var{phi}, @var{g}, @var{h})
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
23 ## @deftypefnx {Function File} {[@dots{}] =} sqp (@var{x0}, @var{phi}, @var{g}, @var{h}, @var{lb}, @var{ub})
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
24 ## @deftypefnx {Function File} {[@dots{}] =} sqp (@var{x0}, @var{phi}, @var{g}, @var{h}, @var{lb}, @var{ub}, @var{maxiter})
13207
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
25 ## @deftypefnx {Function File} {[@dots{}] =} sqp (@var{x0}, @var{phi}, @var{g}, @var{h}, @var{lb}, @var{ub}, @var{maxiter}, @var{tol})
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
26 ## Solve the nonlinear program
6741
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
27 ## @tex
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
28 ## $$
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
29 ## \min_x \phi (x)
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
30 ## $$
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
31 ## @end tex
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
32 ## @ifnottex
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
33 ##
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
34 ## @example
9051
1bf0ce0930be Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents: 8920
diff changeset
35 ## @group
14327
4d917a6a858b doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
36 ## min phi (x)
4d917a6a858b doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
37 ## x
9051
1bf0ce0930be Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents: 8920
diff changeset
38 ## @end group
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
39 ## @end example
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
40 ##
6741
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
41 ## @end ifnottex
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
42 ## subject to
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
43 ## @tex
6741
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
44 ## $$
8167
17352ccd860e describe additional arguments in sqp() documentation string
Ivan Sutoris <ivan.sutoris@gmail.com>
parents: 8047
diff changeset
45 ## 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
46 ## $$
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
47 ## @end tex
6741
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
48 ## @ifnottex
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
49 ##
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
50 ## @example
9051
1bf0ce0930be Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents: 8920
diff changeset
51 ## @group
14327
4d917a6a858b doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
52 ## g(x) = 0
4d917a6a858b doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
53 ## h(x) >= 0
4d917a6a858b doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
54 ## lb <= x <= ub
9051
1bf0ce0930be Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents: 8920
diff changeset
55 ## @end group
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
56 ## @end example
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10678
diff changeset
57 ##
6741
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
58 ## @end ifnottex
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
59 ## @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
60 ## using a sequential quadratic programming method.
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
61 ##
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
62 ## 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
63 ##
7001
8b0cfeb06365 [project @ 2007-10-10 18:02:59 by jwe]
jwe
parents: 6921
diff changeset
64 ## The second argument is a function handle pointing to the objective
13207
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
65 ## function @var{phi}. The objective function must accept one vector
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
66 ## argument and return a scalar.
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
67 ##
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
68 ## The second argument may also be a 2- or 3-element cell array of
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
69 ## function handles. The first element should point to the objective
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
70 ## function, the second should point to a function that computes the
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
71 ## gradient of the objective function, and the third should point to a
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
72 ## function that computes the Hessian of the objective function. If the
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
73 ## gradient function is not supplied, the gradient is computed by finite
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
74 ## differences. If the Hessian function is not supplied, a BFGS update
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
75 ## formula is used to approximate the Hessian.
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
76 ##
13207
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
77 ## When supplied, the gradient function @code{@var{phi}@{2@}} must accept
13931
9de488c6c59c doc: Spellcheck documentation before 3.6.0 release
Rik <octave@nomad.inbox5.com>
parents: 13929
diff changeset
78 ## one vector argument and return a vector. When supplied, the Hessian
13207
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
79 ## function @code{@var{phi}@{3@}} must accept one vector argument and
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
80 ## return a matrix.
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
81 ##
13207
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
82 ## The third and fourth arguments @var{g} and @var{h} are function
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
83 ## handles pointing to functions that compute the equality constraints
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
84 ## and the inequality constraints, respectively. If the problem does
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
85 ## not have equality (or inequality) constraints, then use an empty
13929
9cae456085c2 Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents: 13305
diff changeset
86 ## matrix ([]) for @var{g} (or @var{h}). When supplied, these equality
13207
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
87 ## and inequality constraint functions must accept one vector argument
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
88 ## and return a vector.
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
89 ##
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
90 ## 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
91 ## function handles. The first element should point to the constraint
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
92 ## 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
93 ## gradient of the constraint function:
6741
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
94 ## @tex
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
95 ## $$
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
96 ## \Bigg( {\partial f(x) \over \partial x_1},
6741
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
97 ## {\partial f(x) \over \partial x_2}, \ldots,
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
98 ## {\partial f(x) \over \partial x_N} \Bigg)^T
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
99 ## $$
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
100 ## @end tex
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
101 ## @ifnottex
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10678
diff changeset
102 ##
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
103 ## @example
9051
1bf0ce0930be Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents: 8920
diff changeset
104 ## @group
14327
4d917a6a858b doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
105 ## [ 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
106 ## transpose ( [ ------ ----- ... ------ ] )
4d917a6a858b doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
107 ## [ dx_1 dx_2 dx_N ]
9051
1bf0ce0930be Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents: 8920
diff changeset
108 ## @end group
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
109 ## @end example
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10678
diff changeset
110 ##
6741
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6527
diff changeset
111 ## @end ifnottex
13207
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
112 ## The fifth and sixth arguments, @var{lb} and @var{ub}, contain lower
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
113 ## and upper bounds on @var{x}. These must be consistent with the
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
114 ## equality and inequality constraints @var{g} and @var{h}. If the
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
115 ## arguments are vectors then @var{x}(i) is bound by @var{lb}(i) and
13929
9cae456085c2 Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents: 13305
diff changeset
116 ## @var{ub}(i). A bound can also be a scalar in which case all elements
13207
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
117 ## of @var{x} will share the same bound. If only one bound (lb, ub) is
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
118 ## specified then the other will default to (-@var{realmax},
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
119 ## +@var{realmax}).
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
120 ##
13207
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
121 ## The seventh argument @var{maxiter} specifies the maximum number of
13929
9cae456085c2 Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents: 13305
diff changeset
122 ## iterations. The default value is 100.
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
123 ##
13207
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
124 ## The eighth argument @var{tol} specifies the tolerance for the
14868
5d3a684236b0 maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents: 14595
diff changeset
125 ## stopping criteria. The default value is @code{sqrt (eps)}.
8167
17352ccd860e describe additional arguments in sqp() documentation string
Ivan Sutoris <ivan.sutoris@gmail.com>
parents: 8047
diff changeset
126 ##
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
127 ## The value returned in @var{info} may be one of the following:
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10678
diff changeset
128 ##
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
129 ## @table @asis
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
130 ## @item 101
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
131 ## The algorithm terminated normally.
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
132 ## Either all constraints meet the requested tolerance, or the stepsize,
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
133 ## @tex
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
134 ## $\Delta x,$
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
135 ## @end tex
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
136 ## @ifnottex
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
137 ## delta @var{x},
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
138 ## @end ifnottex
13207
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
139 ## is less than @code{@var{tol} * norm (x)}.
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
140 ##
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
141 ## @item 102
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
142 ## The BFGS update failed.
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
143 ##
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
144 ## @item 103
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
145 ## The maximum number of iterations was reached.
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
146 ## @end table
8167
17352ccd860e describe additional arguments in sqp() documentation string
Ivan Sutoris <ivan.sutoris@gmail.com>
parents: 8047
diff changeset
147 ##
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
148 ## An example of calling @code{sqp}:
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
149 ##
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
150 ## @example
7031
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
151 ## function r = g (x)
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
152 ## r = [ sumsq(x)-10;
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
153 ## 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
154 ## x(1)^3+x(2)^3+1 ];
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
155 ## endfunction
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
156 ##
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
157 ## 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
158 ## 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
159 ## endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
160 ##
7031
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
161 ## 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
162 ##
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
163 ## [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
164 ##
7031
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
165 ## x =
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
166 ##
7031
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
167 ## -1.71714
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
168 ## 1.59571
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
169 ## 1.82725
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
170 ## -0.76364
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
171 ## -0.76364
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
172 ##
7031
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
173 ## obj = 0.053950
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
174 ## info = 101
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
175 ## iter = 8
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
176 ## nf = 10
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
177 ## lambda =
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
178 ##
7031
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
179 ## -0.0401627
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
180 ## 0.0379578
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
181 ## -0.0052227
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
182 ## @end example
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
183 ##
5642
2618a0750ae6 [project @ 2006-03-06 21:26:48 by jwe]
jwe
parents: 5307
diff changeset
184 ## @seealso{qp}
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
185 ## @end deftypefn
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
186
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
187 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
188
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
189 globals = struct (); # data and handles, needed and changed by
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
190 # subfunctions
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
191
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
192 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
193 print_usage ();
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
194 endif
5289
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 if (!isvector (x0))
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
197 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
198 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
199 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
200 x0 = x0';
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
201 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
202
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
203 have_hess = 0;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
204 if (iscell (objf))
11120
a44f979a35ce style fixes for some .m files
John W. Eaton <jwe@octave.org>
parents: 10821
diff changeset
205 switch (numel (objf))
a44f979a35ce style fixes for some .m files
John W. Eaton <jwe@octave.org>
parents: 10821
diff changeset
206 case 1
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
207 obj_fun = objf{1};
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
208 obj_grd = @ (x) fd_obj_grd (x, obj_fun);
11120
a44f979a35ce style fixes for some .m files
John W. Eaton <jwe@octave.org>
parents: 10821
diff changeset
209 case 2
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
210 obj_fun = objf{1};
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
211 obj_grd = objf{2};
11120
a44f979a35ce style fixes for some .m files
John W. Eaton <jwe@octave.org>
parents: 10821
diff changeset
212 case 3
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
213 obj_fun = objf{1};
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
214 obj_grd = objf{2};
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
215 obj_hess = objf{3};
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
216 have_hess = 1;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
217 otherwise
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
218 error ("sqp: invalid objective function specification");
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
219 endswitch
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
220 else
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
221 obj_fun = objf; # No cell array, only obj_fun set
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
222 obj_grd = @ (x) fd_obj_grd (x, obj_fun);
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
223 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
224
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
225 ce_fun = @empty_cf;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
226 ce_grd = @empty_jac;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
227 if (nargin > 2)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
228 if (iscell (cef))
11120
a44f979a35ce style fixes for some .m files
John W. Eaton <jwe@octave.org>
parents: 10821
diff changeset
229 switch (numel (cef))
a44f979a35ce style fixes for some .m files
John W. Eaton <jwe@octave.org>
parents: 10821
diff changeset
230 case 1
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
231 ce_fun = cef{1};
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
232 ce_grd = @ (x) fd_ce_jac (x, ce_fun);
11120
a44f979a35ce style fixes for some .m files
John W. Eaton <jwe@octave.org>
parents: 10821
diff changeset
233 case 2
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
234 ce_fun = cef{1};
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
235 ce_grd = cef{2};
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
236 otherwise
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
237 error ("sqp: invalid equality constraint function specification");
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
238 endswitch
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
239 elseif (! isempty (cef))
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
240 ce_fun = cef; # No cell array, only constraint equality function set
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
241 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
242 endif
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
243 endif
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
244
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
245 ci_fun = @empty_cf;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
246 ci_grd = @empty_jac;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
247 if (nargin > 3)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
248 ## 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
249 globals.cif = cif;
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
250 ## 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
251 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
252 if (iscell (cif))
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
253 if (length (cif) > 0)
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
254 globals.cifcn = cif{1};
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
255 endif
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
256 elseif (! isempty (cif))
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
257 globals.cifcn = cif;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
258 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
259
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
260 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
261 ## constraint inequality function only without any bounds
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
262 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
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 switch length (cif)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
265 case {1}
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
266 ci_fun = cif{1};
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
267 case {2}
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
268 ci_fun = cif{1};
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
269 ci_grd = cif{2};
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
270 otherwise
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
271 error ("sqp: invalid inequality constraint function specification");
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
272 endswitch
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
273 elseif (! isempty (cif))
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
274 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
275 endif
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
276 else
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
277 ## 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
278 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
279 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
280 if (isvector (lb))
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
281 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
282 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
283 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
284 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
285 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
286 if (isa (x0, "single"))
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
287 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
288 else
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
289 globals.lb = tmp_lb = -realmax;
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10542
diff changeset
290 endif
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
291 else
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
292 error ("sqp: invalid lower bound");
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
293 endif
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
294
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
295 if (isvector (ub))
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
296 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
297 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
298 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
299 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
300 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
301 if (isa (x0, "single"))
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
302 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
303 else
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
304 globals.ub = tmp_ub = realmax;
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10542
diff changeset
305 endif
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
306 else
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
307 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
308 endif
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
309
11347
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
310 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
311 error ("sqp: upper bound smaller than lower bound");
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
312 endif
11347
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
313 bounds_grad = [lb_grad; ub_grad];
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
314 ci_fun = @ (x) cf_ub_lb (x, lb_idx, ub_idx, globals);
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
315 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
316 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
317
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
318 endif # if (nargin > 3)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
319
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
320 iter_max = 100;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
321 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
322 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
323 iter_max = maxiter;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
324 else
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
325 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
326 endif
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
327 endif
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
328
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
329 tol = sqrt (eps);
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
330 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
331 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
332 tol = tolerance;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
333 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
334 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
335 endif
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
336 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
337
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
338 ## 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
339 ## 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
340 ## 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
341 ##
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
342 ## obj_fun -- objective function
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
343 ## obj_grad -- objective gradient
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
344 ## 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
345 ## 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
346 ## 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
347 x = x0;
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
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 obj = feval (obj_fun, x0);
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
350 globals.nfun = 1;
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
351
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
352 c = feval (obj_grd, x0);
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
353
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
354 ## 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
355 n = length (x0);
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
356 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
357 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
358 else
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
359 B = eye (n, n);
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
360 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
361
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
362 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
363 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
364
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
365 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
366 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
367
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
368 A = [F; C];
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
369
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
370 ## 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
371 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
372
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
373 qp_iter = 1;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
374 alpha = 1;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
375
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
376 info = 0;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
377 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
378 # 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
379 # report (iter, qp_iter, alpha, __sqp_nfun__, obj);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
380
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
381 while (++iter < iter_max)
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 ## 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
384 ## order necessary conditions.
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
385 nr_f = rows (F);
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
386
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
387 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
388 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
389
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
390 con = [ce; ci];
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
391
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
392 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
393 t1 = norm (ce);
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
394 t2 = all (ci >= 0);
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
395 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
396 t4 = norm (lambda .* con);
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 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
399 info = 101;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
400 break;
6382
f74e71ef6612 [project @ 2007-03-05 21:29:40 by jwe]
jwe
parents: 6046
diff changeset
401 endif
f74e71ef6612 [project @ 2007-03-05 21:29:40 by jwe]
jwe
parents: 6046
diff changeset
402
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
403 ## 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
404 g = -ce;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
405 d = -ci;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
406
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
407 [p, obj_qp, INFO, lambda] = qp (x, B, c, F, g, [], [], d, C,
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
408 Inf (size (d)));
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
409
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
410 info = INFO.info;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
411
13206
658aa9fd8665 sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents: 11587
diff changeset
412 ## FIXME -- check QP solution and attempt to recover if it has
658aa9fd8665 sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents: 11587
diff changeset
413 ## failed. For now, just warn about possible problems.
658aa9fd8665 sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents: 11587
diff changeset
414
658aa9fd8665 sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents: 11587
diff changeset
415 id = "Octave:SQP-QP-subproblem";
658aa9fd8665 sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents: 11587
diff changeset
416 switch (info)
658aa9fd8665 sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents: 11587
diff changeset
417 case 2
658aa9fd8665 sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents: 11587
diff changeset
418 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
419 case 3
658aa9fd8665 sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents: 11587
diff changeset
420 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
421 INFO.solveiter);
658aa9fd8665 sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents: 11587
diff changeset
422 case 6
658aa9fd8665 sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents: 11587
diff changeset
423 warning (id, "sqp: QP subproblem is infeasible");
658aa9fd8665 sqp: provide warnings for QP subproblem failures
John W. Eaton <jwe@octave.org>
parents: 11587
diff changeset
424 endswitch
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
425
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
426 ## 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
427 ## merit function phi.
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
428 [x_new, alpha, obj_new, globals] = \
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
429 linesearch_L1 (x, p, obj_fun, obj_grd, ce_fun, ci_fun, lambda, \
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
430 obj, globals);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
431
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
432 ## Evaluate objective function, constraints, and gradients at x_new.
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
433 c_new = feval (obj_grd, x_new);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
434
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
435 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
436 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
437
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
438 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
439 C_new = feval (ci_grd, x_new);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
440
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
441 A_new = [F_new; C_new];
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
442
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
443 ## Set
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 ## s = alpha * p
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
446 ## 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
447
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
448 y = c_new - c;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
449
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
450 if (! isempty (A))
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
451 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
452 y -= t;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
453 endif
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
454
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
455 delx = x_new - x;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
456
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
457 if (norm (delx) < tol * norm (x))
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
458 info = 101;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
459 break;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
460 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
461
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
462 if (have_hess)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
463
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
464 B = feval (obj_hess, x);
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 else
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
467 ## 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
468 delxt = delx';
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
469
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
470 ## 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
471 ## 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
472 d1 = delxt*B*delx;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
473
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
474 t1 = 0.2 * d1;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
475 t2 = delxt*y;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
476
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
477 if (t2 < t1)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
478 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
479 else
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
480 theta = 1;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
481 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
482
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
483 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
484
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
485 d2 = delxt*r;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
486
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
487 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
488 info = 102;
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10542
diff changeset
489 break;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
490 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
491
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
492 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
493
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
494 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
495
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
496 x = x_new;
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 obj = obj_new;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
499
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
500 c = c_new;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
501
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
502 ce = ce_new;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
503 F = F_new;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
504
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
505 ci = ci_new;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
506 C = C_new;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
507
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
508 A = A_new;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
509
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
510 # 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
511
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
512 endwhile
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
513
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
514 if (iter >= iter_max)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
515 info = 103;
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
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
518 nf = globals.nfun;
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
519
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
520 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
521
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
522
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
523 function [merit, obj, globals] = phi_L1 (obj, obj_fun, ce_fun, ci_fun, \
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
524 x, mu, globals)
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
525
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
526 ce = feval (ce_fun, x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
527 ci = feval (ci_fun, x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
528
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
529 idx = ci < 0;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
530
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
531 con = [ce; ci(idx)];
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
532
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
533 if (isempty (obj))
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
534 obj = feval (obj_fun, x);
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
535 globals.nfun++;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
536 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
537
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
538 merit = obj;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
539 t = norm (con, 1) / mu;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
540
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
541 if (! isempty (t))
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
542 merit += t;
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
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
545 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
546
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
547
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
548 function [x_new, alpha, obj, globals] = \
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
549 linesearch_L1 (x, p, obj_fun, obj_grd, ce_fun, ci_fun, lambda, \
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
550 obj, globals)
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
551
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
552 ## Choose parameters
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
553 ##
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
554 ## eta in the range (0, 0.5)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
555 ## tau in the range (0, 1)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
556
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
557 eta = 0.25;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
558 tau = 0.5;
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 delta_bar = sqrt (eps);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
561
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
562 if (isempty (lambda))
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
563 mu = 1 / delta_bar;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
564 else
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
565 mu = 1 / (norm (lambda, Inf) + delta_bar);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
566 endif
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 alpha = 1;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
569
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
570 c = feval (obj_grd, x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
571 ce = feval (ce_fun, x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
572
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
573 [phi_x_mu, obj, globals] = phi_L1 (obj, obj_fun, ce_fun, ci_fun, x, \
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
574 mu, globals);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
575
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
576 D_phi_x_mu = c' * p;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
577 d = feval (ci_fun, x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
578 ## only those elements of d corresponding
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
579 ## to violated constraints should be included.
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
580 idx = d < 0;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
581 t = - norm ([ce; d(idx)], 1) / mu;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
582 if (! isempty (t))
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
583 D_phi_x_mu += t;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
584 endif
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 while (1)
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
587 [p1, obj, globals] = phi_L1 ([], obj_fun, ce_fun, ci_fun, \
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
588 x+alpha*p, mu, globals);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
589 p2 = phi_x_mu+eta*alpha*D_phi_x_mu;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
590 if (p1 > p2)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
591 ## 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
592 ## range (0, tau).
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
593 tau_alpha = 0.9 * tau; # ??
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
594 alpha = tau_alpha * alpha;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
595 else
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
596 break;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
597 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
598 endwhile
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
599
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
600 x_new = x + alpha * p;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
601
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
602 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
603
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
604
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
605 function grd = fdgrd (f, x)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
606
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
607 if (! isempty (f))
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
608 y0 = feval (f, x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
609 nx = length (x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
610 grd = zeros (nx, 1);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
611 deltax = sqrt (eps);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
612 for i = 1:nx
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
613 t = x(i);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
614 x(i) += deltax;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
615 grd(i) = (feval (f, x) - y0) / deltax;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
616 x(i) = t;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
617 endfor
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
618 else
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
619 grd = zeros (0, 1);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
620 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
621
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
622 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
623
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 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
626
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
627 nx = length (x);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
628 if (! isempty (f))
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
629 y0 = feval (f, x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
630 nf = length (y0);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
631 nx = length (x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
632 jac = zeros (nf, nx);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
633 deltax = sqrt (eps);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
634 for i = 1:nx
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
635 t = x(i);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
636 x(i) += deltax;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
637 jac(:,i) = (feval (f, x) - y0) / deltax;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
638 x(i) = t;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
639 endfor
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
640 else
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
641 jac = zeros (0, nx);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
642 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
643
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
644 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
645
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
646
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
647 function grd = fd_obj_grd (x, obj_fun)
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
648
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
649 grd = fdgrd (obj_fun, x);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
650
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
651 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
652
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
653
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
654 function res = empty_cf (x)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
655
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
656 res = zeros (0, 1);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
657
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
658 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
659
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
660
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
661 function res = empty_jac (x)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
662
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
663 res = zeros (0, length (x));
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
664
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
665 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
666
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
667
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
668 function jac = fd_ce_jac (x, ce_fun)
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
669
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
670 jac = fdjac (ce_fun, x);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
671
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
672 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
673
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
674
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
675 function jac = fd_ci_jac (x, cifcn)
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
676
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
677 ## 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
678 jac = fdjac (cifcn, x);
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
679
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
680 endfunction
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
681
7017
a1dbe9d80eee [project @ 2007-10-12 21:27:11 by jwe]
jwe
parents: 7016
diff changeset
682
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
683 function res = cf_ub_lb (x, lbidx, ubidx, globals)
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
684
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
685 ## 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
686 if (isempty (globals.cifcn))
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
687 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
688 else
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
689 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
690 globals.ub-x(ubidx,1)];
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
691 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
692
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
693 endfunction
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
694
7017
a1dbe9d80eee [project @ 2007-10-12 21:27:11 by jwe]
jwe
parents: 7016
diff changeset
695
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
696 function res = cigrad_ub_lb (x, bgrad, globals)
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
697
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
698 cigradfcn = @ (x) fd_ci_jac (x, globals.cifcn);
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
699
14595
11a9d448fdc3 Get rid of global variables in sqp.
Olaf Till <i7tiol@t-online.de>
parents: 14552
diff changeset
700 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
701 cigradfcn = globals.cif{2};
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
702 endif
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
703
11347
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
704 if (isempty (cigradfcn))
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
705 res = bgrad;
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
706 else
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
707 res = [feval(cigradfcn,x); bgrad];
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
708 endif
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
709
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
710 endfunction
7361
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
711
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
712 # Utility function used to debug sqp
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
713 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
714
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
715 if (nargin == 0)
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
716 printf (" Itn ItQP Step Nfun Objective\n");
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
717 else
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
718 printf ("%5d %4d %8.1g %5d %13.6e\n", 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
719 endif
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
720
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
721 endfunction
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
722
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14327
diff changeset
723
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
724 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
725 %% Test Code
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
726
13305
63463570d9fe Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents: 13268
diff changeset
727 %!function r = __g (x)
7371
fe9a44d753d6 [project @ 2008-01-14 19:34:22 by jwe]
jwe
parents: 7361
diff changeset
728 %! r = [sumsq(x)-10;
fe9a44d753d6 [project @ 2008-01-14 19:34:22 by jwe]
jwe
parents: 7361
diff changeset
729 %! x(2)*x(3)-5*x(4)*x(5);
fe9a44d753d6 [project @ 2008-01-14 19:34:22 by jwe]
jwe
parents: 7361
diff changeset
730 %! 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
731 %!endfunction
7361
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
732 %!
13305
63463570d9fe Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents: 13268
diff changeset
733 %!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
734 %! 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
735 %!endfunction
7361
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
736 %!
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
737 %!test
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
738 %!
7361
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
739 %! 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
740 %!
13305
63463570d9fe Add %!endfunction block keyword to test.m
Rik <octave@nomad.inbox5.com>
parents: 13268
diff changeset
741 %! [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
742 %!
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
743 %! x_opt = [-1.717143501952599;
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
744 %! 1.595709610928535;
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
745 %! 1.827245880097156;
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
746 %! -0.763643103133572;
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
747 %! -0.763643068453300];
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
748 %!
7371
fe9a44d753d6 [project @ 2008-01-14 19:34:22 by jwe]
jwe
parents: 7361
diff changeset
749 %! obj_opt = 0.0539498477702739;
7361
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
750 %!
15757
534a2c881f45 Increase test tolerance on sqp (bug #37742)
Mike Miller <mtmiller@ieee.org>
parents: 14868
diff changeset
751 %! 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
752 %! 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
753
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
754 %% Test input validation
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
755 %!error sqp ()
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
756 %!error sqp (1)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
757 %!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
758 %!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
759 %!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
760 %!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
761 %!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
762 %!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
763 %!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
764 %!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
765 %!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
766 %!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
767 %!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
768 %!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
769 %!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
770 %!error sqp (1, cell (3,1), cell (2,1), cell (2,1),[],[],[],-1)