annotate scripts/optimization/sqp.m @ 13207:88bd1d1d6657 stable

Reword sqp's docstring
author Jordi Gutiérrez Hermoso <jordigh@octave.org>
date Sat, 24 Sep 2011 11:04:24 -0500
parents c792872f8942
children 571557ddabb9
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
11523
fd0a3ac60b0e update copyright notices
John W. Eaton <jwe@octave.org>
parents: 11472
diff changeset
1 ## Copyright (C) 2005-2011 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
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
36 ## min phi (x)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
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
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
52 ## g(x) = 0
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
53 ## h(x) >= 0
8167
17352ccd860e describe additional arguments in sqp() documentation string
Ivan Sutoris <ivan.sutoris@gmail.com>
parents: 8047
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
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
78 ## one vector argument and return a vector. When supplifed, the Hessian
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
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
86 ## matrix ([]) for @var{g} (or @var{h}). When supplied, these equality
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
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
105 ## [ d f(x) d f(x) d f(x) ]
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
106 ## transpose ( [ ------ ----- ... ------ ] )
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
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
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
diff changeset
116 ## @var{ub}(i). A bound can also be a scalar in which case all elements
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
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
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
88bd1d1d6657 Reword sqp's docstring
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 11587
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)
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
158 ## obj = exp(prod(x)) - 0.5*(x(1)^3+x(2)^3+1)^2;
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
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
189 global __sqp_nfun__;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
190 global __sqp_obj_fun__;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
191 global __sqp_ce_fun__;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
192 global __sqp_ci_fun__;
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
193 global __sqp_cif__;
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
194 global __sqp_cifcn__;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
195
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
196 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
197 print_usage ();
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
198 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
199
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
200 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
201 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
202 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
203 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
204 x0 = x0';
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
205 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
206
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
207 obj_grd = @fd_obj_grd;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
208 have_hess = 0;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
209 if (iscell (objf))
11120
a44f979a35ce style fixes for some .m files
John W. Eaton <jwe@octave.org>
parents: 10821
diff changeset
210 switch (numel (objf))
a44f979a35ce style fixes for some .m files
John W. Eaton <jwe@octave.org>
parents: 10821
diff changeset
211 case 1
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
212 obj_fun = objf{1};
11120
a44f979a35ce style fixes for some .m files
John W. Eaton <jwe@octave.org>
parents: 10821
diff changeset
213 case 2
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
214 obj_fun = objf{1};
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
215 obj_grd = objf{2};
11120
a44f979a35ce style fixes for some .m files
John W. Eaton <jwe@octave.org>
parents: 10821
diff changeset
216 case 3
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
217 obj_fun = objf{1};
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
218 obj_grd = objf{2};
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
219 obj_hess = objf{3};
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
220 have_hess = 1;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
221 otherwise
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
222 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
223 endswitch
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
224 else
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
225 obj_fun = objf; # No cell array, only obj_fun set
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
226 endif
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
227 __sqp_obj_fun__ = obj_fun;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
228
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
229 ce_fun = @empty_cf;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
230 ce_grd = @empty_jac;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
231 if (nargin > 2)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
232 ce_grd = @fd_ce_jac;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
233 if (iscell (cef))
11120
a44f979a35ce style fixes for some .m files
John W. Eaton <jwe@octave.org>
parents: 10821
diff changeset
234 switch (numel (cef))
a44f979a35ce style fixes for some .m files
John W. Eaton <jwe@octave.org>
parents: 10821
diff changeset
235 case 1
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
236 ce_fun = cef{1};
11120
a44f979a35ce style fixes for some .m files
John W. Eaton <jwe@octave.org>
parents: 10821
diff changeset
237 case 2
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
238 ce_fun = cef{1};
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
239 ce_grd = cef{2};
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
240 otherwise
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
241 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
242 endswitch
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
243 elseif (! isempty (cef))
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
244 ce_fun = cef; # No cell array, only constraint equality function set
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
245 endif
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
246 endif
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
247 __sqp_ce_fun__ = ce_fun;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
248
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
249 ci_fun = @empty_cf;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
250 ci_grd = @empty_jac;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
251 if (nargin > 3)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
252 ## constraint function given by user with possible gradient
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
253 __sqp_cif__ = cif;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
254 ## constraint function given by user without gradient
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
255 __sqp_cifcn__ = @empty_cf;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
256 if (iscell (cif))
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
257 if (length (cif) > 0)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
258 __sqp_cifcn__ = cif{1};
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
259 endif
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
260 elseif (! isempty (cif))
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
261 __sqp_cifcn__ = cif;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
262 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
263
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
264 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
265 ## constraint inequality function only without any bounds
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
266 ci_grd = @fd_ci_jac;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
267 if (iscell (cif))
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
268 switch length (cif)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
269 case {1}
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
270 ci_fun = cif{1};
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
271 case {2}
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
272 ci_fun = cif{1};
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
273 ci_grd = cif{2};
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
274 otherwise
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
275 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
276 endswitch
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
277 elseif (! isempty (cif))
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
278 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
279 endif
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
280 else
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
281 ## constraint inequality function with bounds present
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
282 global __sqp_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
283 lb_idx = ub_idx = true (size (x0));
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
284 ub_grad = - (lb_grad = eye (rows (x0)));
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
285 if (isvector (lb))
11347
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
286 __sqp_lb__ = tmp_lb = lb(:);
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
287 lb_idx(:) = tmp_idx = (lb != -Inf);
11371
c767bb1afa03 sqp.m: Fix indexing error in sqp bounds selection
Olaf Till <olaf.till@uni-jena.de>
parents: 11348
diff changeset
288 __sqp_lb__ = __sqp_lb__(tmp_idx, 1);
11347
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
289 lb_grad = lb_grad(lb_idx, :);
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
290 elseif (isempty (lb))
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
291 if (isa (x0, "single"))
11347
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
292 __sqp_lb__ = tmp_lb = -realmax ("single");
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
293 else
11347
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
294 __sqp_lb__ = tmp_lb = -realmax;
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10542
diff changeset
295 endif
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
296 else
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
297 error ("sqp: invalid lower bound");
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
298 endif
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
299
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
300 global __sqp_ub__;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
301 if (isvector (ub))
11347
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
302 __sqp_ub__ = tmp_ub = ub(:);
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
303 ub_idx(:) = tmp_idx = (ub != Inf);
11371
c767bb1afa03 sqp.m: Fix indexing error in sqp bounds selection
Olaf Till <olaf.till@uni-jena.de>
parents: 11348
diff changeset
304 __sqp_ub__ = __sqp_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
305 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
306 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
307 if (isa (x0, "single"))
11347
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
308 __sqp_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
309 else
11347
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
310 __sqp_ub__ = tmp_ub = realmax;
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10542
diff changeset
311 endif
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
312 else
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
313 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
314 endif
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
315
11347
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
316 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
317 error ("sqp: upper bound smaller than lower bound");
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
318 endif
11347
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
319 bounds_grad = [lb_grad; ub_grad];
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
320 ci_fun = @ (x) cf_ub_lb (x, lb_idx, ub_idx);
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
321 ci_grd = @ (x) cigrad_ub_lb (x, bounds_grad);
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
322 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
323
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
324 __sqp_ci_fun__ = ci_fun;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
325 endif # if (nargin > 3)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
326
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
327 iter_max = 100;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
328 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
329 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
330 iter_max = maxiter;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
331 else
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
332 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
333 endif
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
334 endif
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
335
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
336 tol = sqrt (eps);
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
337 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
338 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
339 tol = tolerance;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
340 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
341 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
342 endif
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
343 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
344
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
345 ## 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
346 ## 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
347 ## 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
348 ##
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
349 ## obj_fun -- objective function
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
350 ## obj_grad -- objective gradient
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
351 ## 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
352 ## 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
353 ## 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
354 x = x0;
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
355
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
356 obj = feval (obj_fun, x0);
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
357 __sqp_nfun__ = 1;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
358
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
359 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
360
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
361 ## 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
362 n = length (x0);
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
363 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
364 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
365 else
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
366 B = eye (n, n);
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
367 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
368
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
369 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
370 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
371
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
372 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
373 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
374
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
375 A = [F; C];
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
376
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
377 ## 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
378 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
379
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
380 qp_iter = 1;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
381 alpha = 1;
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 info = 0;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
384 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
385 # 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
386 # report (iter, qp_iter, alpha, __sqp_nfun__, obj);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
387
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
388 while (++iter < iter_max)
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 ## 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
391 ## order necessary conditions.
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
392 nr_f = rows (F);
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
393
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
394 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
395 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
396
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
397 con = [ce; ci];
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
398
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
399 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
400 t1 = norm (ce);
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
401 t2 = all (ci >= 0);
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
402 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
403 t4 = norm (lambda .* con);
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
404
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
405 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
406 info = 101;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
407 break;
6382
f74e71ef6612 [project @ 2007-03-05 21:29:40 by jwe]
jwe
parents: 6046
diff changeset
408 endif
f74e71ef6612 [project @ 2007-03-05 21:29:40 by jwe]
jwe
parents: 6046
diff changeset
409
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
410 ## 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
411 g = -ce;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
412 d = -ci;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
413
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
414 [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
415 Inf (size (d)));
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
416
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
417 info = INFO.info;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
418
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
419 ## Check QP solution and attempt to recover if it has failed.
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
420
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
421 ## 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
422 ## merit function phi.
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
423 [x_new, alpha, obj_new] = linesearch_L1 (x, p, obj_fun, obj_grd,
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
424 ce_fun, ci_fun, lambda, obj);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
425
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
426 ## 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
427 c_new = feval (obj_grd, x_new);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
428
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
429 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
430 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
431
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
432 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
433 C_new = feval (ci_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 A_new = [F_new; C_new];
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
436
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
437 ## Set
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
438 ##
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
439 ## s = alpha * p
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
440 ## 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
441
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
442 y = c_new - c;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
443
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
444 if (! isempty (A))
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
445 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
446 y -= t;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
447 endif
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
448
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
449 delx = x_new - x;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
450
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
451 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
452 info = 101;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
453 break;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
454 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
455
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
456 if (have_hess)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
457
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
458 B = feval (obj_hess, x);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
459
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
460 else
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
461 ## 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
462 delxt = delx';
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 ## 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
465 ## 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
466 d1 = delxt*B*delx;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
467
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
468 t1 = 0.2 * d1;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
469 t2 = delxt*y;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
470
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
471 if (t2 < t1)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
472 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
473 else
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
474 theta = 1;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
475 endif
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 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
478
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
479 d2 = delxt*r;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
480
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
481 if (d1 == 0 || d2 == 0)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
482 info = 102;
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10542
diff changeset
483 break;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
484 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
485
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
486 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
487
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
488 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
489
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
490 x = x_new;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
491
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
492 obj = obj_new;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
493
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
494 c = c_new;
5289
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 ce = ce_new;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
497 F = F_new;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
498
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
499 ci = ci_new;
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 A = A_new;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
503
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
504 # 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
505
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
506 endwhile
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 if (iter >= iter_max)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
509 info = 103;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
510 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
511
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
512 nf = __sqp_nfun__;
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
513
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
514 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
515
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
516
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
517 function [merit, obj] = phi_L1 (obj, obj_fun, ce_fun, ci_fun, x, mu)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
518
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
519 global __sqp_nfun__;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
520
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
521 ce = feval (ce_fun, x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
522 ci = feval (ci_fun, x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
523
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
524 idx = ci < 0;
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 con = [ce; ci(idx)];
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
527
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
528 if (isempty (obj))
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
529 obj = feval (obj_fun, x);
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
530 __sqp_nfun__++;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
531 endif
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 merit = obj;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
534 t = norm (con, 1) / mu;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
535
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
536 if (! isempty (t))
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
537 merit += t;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
538 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
539
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
540 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
541
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
542
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
543 function [x_new, alpha, obj] = linesearch_L1 (x, p, obj_fun, obj_grd,
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10542
diff changeset
544 ce_fun, ci_fun, lambda, obj)
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
545
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
546 ## Choose parameters
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
547 ##
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
548 ## eta in the range (0, 0.5)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
549 ## tau in the range (0, 1)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
550
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
551 eta = 0.25;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
552 tau = 0.5;
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 delta_bar = sqrt (eps);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
555
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
556 if (isempty (lambda))
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
557 mu = 1 / delta_bar;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
558 else
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
559 mu = 1 / (norm (lambda, Inf) + delta_bar);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
560 endif
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 alpha = 1;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
563
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
564 c = feval (obj_grd, x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
565 ce = feval (ce_fun, x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
566
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
567 [phi_x_mu, obj] = phi_L1 (obj, obj_fun, ce_fun, ci_fun, x, mu);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
568
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
569 D_phi_x_mu = c' * p;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
570 d = feval (ci_fun, x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
571 ## only those elements of d corresponding
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
572 ## to violated constraints should be included.
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
573 idx = d < 0;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
574 t = - norm ([ce; d(idx)], 1) / mu;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
575 if (! isempty (t))
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
576 D_phi_x_mu += t;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
577 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
578
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
579 while (1)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
580 [p1, obj] = phi_L1 ([], obj_fun, ce_fun, ci_fun, x+alpha*p, mu);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
581 p2 = phi_x_mu+eta*alpha*D_phi_x_mu;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
582 if (p1 > p2)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
583 ## 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
584 ## range (0, tau).
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
585 tau_alpha = 0.9 * tau; # ??
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
586 alpha = tau_alpha * alpha;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
587 else
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
588 break;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
589 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
590 endwhile
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
591
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
592 x_new = x + alpha * p;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
593
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
594 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
595
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
596
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
597 function grd = fdgrd (f, x)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
598
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
599 if (! isempty (f))
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
600 y0 = feval (f, x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
601 nx = length (x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
602 grd = zeros (nx, 1);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
603 deltax = sqrt (eps);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
604 for i = 1:nx
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
605 t = x(i);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
606 x(i) += deltax;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
607 grd(i) = (feval (f, x) - y0) / deltax;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
608 x(i) = t;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
609 endfor
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
610 else
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
611 grd = zeros (0, 1);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
612 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
613
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
614 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
615
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
616
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
617 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
618
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
619 nx = length (x);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
620 if (! isempty (f))
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
621 y0 = feval (f, x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
622 nf = length (y0);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
623 nx = length (x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
624 jac = zeros (nf, nx);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
625 deltax = sqrt (eps);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
626 for i = 1:nx
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
627 t = x(i);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
628 x(i) += deltax;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
629 jac(:,i) = (feval (f, x) - y0) / deltax;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
630 x(i) = t;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
631 endfor
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
632 else
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
633 jac = zeros (0, nx);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
634 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
635
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
636 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
637
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
638
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
639 function grd = fd_obj_grd (x)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
640
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
641 global __sqp_obj_fun__;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
642
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
643 grd = fdgrd (__sqp_obj_fun__, x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
644
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
645 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
646
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
647
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
648 function res = empty_cf (x)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
649
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
650 res = zeros (0, 1);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
651
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
652 endfunction
5289
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
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
655 function res = empty_jac (x)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
656
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
657 res = zeros (0, length (x));
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
658
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
659 endfunction
5289
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
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
662 function jac = fd_ce_jac (x)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
663
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
664 global __sqp_ce_fun__;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
665
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
666 jac = fdjac (__sqp_ce_fun__, x);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
667
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
668 endfunction
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
669
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
670
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
671 function jac = fd_ci_jac (x)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
672
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
673 global __sqp_cifcn__;
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
674 ## __sqp_cifcn__ = constraint function without gradients and lb or ub
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
675 jac = fdjac (__sqp_cifcn__, x);
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
676
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
677 endfunction
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
678
7017
a1dbe9d80eee [project @ 2007-10-12 21:27:11 by jwe]
jwe
parents: 7016
diff changeset
679
11347
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
680 function res = cf_ub_lb (x, lbidx, ubidx)
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
681
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
682 ## combine constraint function with ub and lb
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
683 global __sqp_cifcn__ __sqp_lb__ __sqp_ub__
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
684
11347
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
685 if (isempty (__sqp_cifcn__))
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
686 res = [x(lbidx,1)-__sqp_lb__; __sqp_ub__-x(ubidx,1)];
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
687 else
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
688 res = [feval(__sqp_cifcn__,x); \
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
689 x(lbidx,1)-__sqp_lb__; __sqp_ub__-x(ubidx,1)];
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
690 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
691
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
692 endfunction
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
693
7017
a1dbe9d80eee [project @ 2007-10-12 21:27:11 by jwe]
jwe
parents: 7016
diff changeset
694
11347
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
695 function res = cigrad_ub_lb (x, bgrad)
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
696
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
697 global __sqp_cif__
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
698
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
699 cigradfcn = @fd_ci_jac;
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
700
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
701 if (iscell (__sqp_cif__) && length (__sqp_cif__) > 1)
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
702 cigradfcn = __sqp_cif__{2};
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
703 endif
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
704
11347
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
705 if (isempty (cigradfcn))
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
706 res = bgrad;
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
707 else
2726132f77f6 sqp.m: Remove never violated Inf bounds from computation (bug #31742)
Rik <octave@nomad.inbox5.com>
parents: 11290
diff changeset
708 res = [feval(cigradfcn,x); bgrad];
6768
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
709 endif
40e1255eda0e [project @ 2007-06-29 20:11:58 by jwe]
jwe
parents: 6741
diff changeset
710
7399
c1a3d6c7d2fb [project @ 2008-01-18 07:37:35 by jwe]
jwe
parents: 7381
diff changeset
711 endfunction
7361
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
712
11348
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
713 # 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
714 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
715
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
716 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
717 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
718 else
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
719 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
720 endif
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
721
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
722 endfunction
2ae0ca4ee36b sqp.m: Change docstring to refer to x0 as the initial seed vector
Rik <octave@nomad.inbox5.com>
parents: 11347
diff changeset
723
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
7371
fe9a44d753d6 [project @ 2008-01-14 19:34:22 by jwe]
jwe
parents: 7361
diff changeset
727 %!function r = g (x)
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 ];
7361
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
731 %!
7371
fe9a44d753d6 [project @ 2008-01-14 19:34:22 by jwe]
jwe
parents: 7361
diff changeset
732 %!function obj = phi (x)
fe9a44d753d6 [project @ 2008-01-14 19:34:22 by jwe]
jwe
parents: 7361
diff changeset
733 %! obj = exp(prod(x)) - 0.5*(x(1)^3+x(2)^3+1)^2;
7361
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
734 %!
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
735 %!test
10678
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
736 %!
7361
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
737 %! 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
738 %!
7381
90b931a40617 [project @ 2008-01-15 19:12:14 by jwe]
jwe
parents: 7371
diff changeset
739 %! [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
740 %!
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
741 %! x_opt = [-1.717143501952599;
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
742 %! 1.595709610928535;
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
743 %! 1.827245880097156;
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
744 %! -0.763643103133572;
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
745 %! -0.763643068453300];
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
746 %!
7371
fe9a44d753d6 [project @ 2008-01-14 19:34:22 by jwe]
jwe
parents: 7361
diff changeset
747 %! obj_opt = 0.0539498477702739;
7361
a2870fd8ac58 [project @ 2008-01-12 07:50:54 by jwe]
jwe
parents: 7031
diff changeset
748 %!
8047
d54f113aa983 Increase test script tolerances.
Thomas Treichl <Thomas.Treichl@gmx.net>
parents: 7795
diff changeset
749 %! assert (all (abs (x-x_opt) < 5*sqrt (eps)) && abs (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
750
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
751 %% Test input validation
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
752 %!error sqp ()
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
753 %!error sqp (1)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
754 %!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
755 %!error sqp (1,2,3,4,5)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
756 %!error sqp (ones(2,2))
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
757 %!error sqp (1,cell(4,1))
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
758 %!error sqp (1,cell(3,1),cell(3,1))
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
759 %!error sqp (1,cell(3,1),cell(2,1),cell(3,1))
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
760 %!error sqp (1,cell(3,1),cell(2,1),cell(2,1),ones(2,2),[])
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
761 %!error sqp (1,cell(3,1),cell(2,1),cell(2,1),[],ones(2,2))
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
762 %!error sqp (1,cell(3,1),cell(2,1),cell(2,1),1,-1)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
763 %!error sqp (1,cell(3,1),cell(2,1),cell(2,1),[],[],ones(2,2))
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
764 %!error sqp (1,cell(3,1),cell(2,1),cell(2,1),[],[],-1)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
765 %!error sqp (1,cell(3,1),cell(2,1),cell(2,1),[],[],1.5)
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
766 %!error sqp (1,cell(3,1),cell(2,1),cell(2,1),[],[],[],ones(2,2))
35338deff753 Guarantee equivalent results if sqp called with or wihout bounds
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
767 %!error sqp (1,cell(3,1),cell(2,1),cell(2,1),[],[],[],-1)