annotate scripts/optimization/qp.m @ 20165:f1d0f506ee78 stable

doc: Update more docstrings to have one sentence summary as first line. Reviewed optimization, polynomial, signal script directories. * scripts/optimization/fminbnd.m, scripts/optimization/fminsearch.m, scripts/optimization/fminunc.m, scripts/optimization/fsolve.m, scripts/optimization/fzero.m, scripts/optimization/glpk.m, scripts/optimization/lsqnonneg.m, scripts/optimization/pqpnonneg.m, scripts/optimization/qp.m, scripts/optimization/sqp.m, scripts/polynomial/compan.m, scripts/polynomial/mkpp.m, scripts/polynomial/mpoles.m, scripts/polynomial/pchip.m, scripts/polynomial/poly.m, scripts/polynomial/polyaffine.m, scripts/polynomial/polyder.m, scripts/polynomial/polyeig.m, scripts/polynomial/polyfit.m, scripts/polynomial/polygcd.m, scripts/polynomial/polyint.m, scripts/polynomial/polyout.m, scripts/polynomial/polyval.m, scripts/polynomial/ppder.m, scripts/polynomial/ppint.m, scripts/polynomial/ppjumps.m, scripts/polynomial/ppval.m, scripts/polynomial/residue.m, scripts/polynomial/roots.m, scripts/polynomial/spline.m, scripts/polynomial/splinefit.m, scripts/polynomial/unmkpp.m, scripts/signal/arch_fit.m, scripts/signal/arch_rnd.m, scripts/signal/arma_rnd.m, scripts/signal/autoreg_matrix.m, scripts/signal/bartlett.m, scripts/signal/blackman.m, scripts/signal/detrend.m, scripts/signal/diffpara.m, scripts/signal/durbinlevinson.m, scripts/signal/fftconv.m, scripts/signal/fftfilt.m, scripts/signal/fftshift.m, scripts/signal/filter2.m, scripts/signal/freqz.m, scripts/signal/hamming.m, scripts/signal/hanning.m, scripts/signal/hurst.m, scripts/signal/ifftshift.m, scripts/signal/periodogram.m, scripts/signal/sinc.m, scripts/signal/sinetone.m, scripts/signal/sinewave.m, scripts/signal/spectral_adf.m, scripts/signal/spectral_xdf.m, scripts/signal/spencer.m, scripts/signal/stft.m, scripts/signal/synthesis.m, scripts/signal/unwrap.m, scripts/signal/yulewalker.m: Update more docstrings to have one sentence summary as first line.
author Rik <rik@octave.org>
date Mon, 04 May 2015 21:50:57 -0700
parents 9fc020886ae9
children
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
19697
4197fc428c7d maint: Update copyright notices for 2015.
John W. Eaton <jwe@octave.org>
parents: 17931
diff changeset
1 ## Copyright (C) 2013-2015 Julien Bect
4197fc428c7d maint: Update copyright notices for 2015.
John W. Eaton <jwe@octave.org>
parents: 17931
diff changeset
2 ## Copyright (C) 2000-2015 Gabriele Pannocchia.
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
3 ##
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
4 ## This file is part of Octave.
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
5 ##
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
6 ## 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
7 ## 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: 6848
diff changeset
8 ## the Free Software Foundation; either version 3 of the License, or (at
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 6848
diff changeset
9 ## your option) any later version.
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
10 ##
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
11 ## 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
12 ## WITHOUT ANY WARRANTY; without even the implied warranty of
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
13 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
14 ## General Public License for more details.
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
15 ##
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
16 ## 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: 6848
diff changeset
17 ## along with Octave; see the file COPYING. If not, see
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 6848
diff changeset
18 ## <http://www.gnu.org/licenses/>.
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
19
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
20 ## -*- texinfo -*-
10793
be55736a0783 Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
21 ## @deftypefn {Function File} {[@var{x}, @var{obj}, @var{info}, @var{lambda}] =} qp (@var{x0}, @var{H})
10059
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
22 ## @deftypefnx {Function File} {[@var{x}, @var{obj}, @var{info}, @var{lambda}] =} qp (@var{x0}, @var{H}, @var{q})
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
23 ## @deftypefnx {Function File} {[@var{x}, @var{obj}, @var{info}, @var{lambda}] =} qp (@var{x0}, @var{H}, @var{q}, @var{A}, @var{b})
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
24 ## @deftypefnx {Function File} {[@var{x}, @var{obj}, @var{info}, @var{lambda}] =} qp (@var{x0}, @var{H}, @var{q}, @var{A}, @var{b}, @var{lb}, @var{ub})
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
25 ## @deftypefnx {Function File} {[@var{x}, @var{obj}, @var{info}, @var{lambda}] =} qp (@var{x0}, @var{H}, @var{q}, @var{A}, @var{b}, @var{lb}, @var{ub}, @var{A_lb}, @var{A_in}, @var{A_ub})
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
26 ## @deftypefnx {Function File} {[@var{x}, @var{obj}, @var{info}, @var{lambda}] =} qp (@dots{}, @var{options})
20165
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
27 ## Solve a quadratic program (QP).
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
28 ##
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
29 ## Solve the quadratic program defined by
6741
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6526
diff changeset
30 ## @tex
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6526
diff changeset
31 ## $$
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6526
diff changeset
32 ## \min_x {1 \over 2} x^T H x + x^T q
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6526
diff changeset
33 ## $$
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6526
diff changeset
34 ## @end tex
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6526
diff changeset
35 ## @ifnottex
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
36 ##
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
37 ## @example
9051
1bf0ce0930be Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents: 8920
diff changeset
38 ## @group
14327
4d917a6a858b doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
39 ## min 0.5 x'*H*x + x'*q
4d917a6a858b doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
40 ## x
9051
1bf0ce0930be Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents: 8920
diff changeset
41 ## @end group
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
42 ## @end example
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
43 ##
6741
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6526
diff changeset
44 ## @end ifnottex
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6526
diff changeset
45 ## subject to
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
46 ## @tex
6741
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6526
diff changeset
47 ## $$
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6526
diff changeset
48 ## Ax = b \qquad lb \leq x \leq ub \qquad A_{lb} \leq A_{in} \leq A_{ub}
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6526
diff changeset
49 ## $$
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
50 ## @end tex
6741
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6526
diff changeset
51 ## @ifnottex
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
52 ##
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
53 ## @example
9051
1bf0ce0930be Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents: 8920
diff changeset
54 ## @group
14327
4d917a6a858b doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
55 ## A*x = b
4d917a6a858b doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
56 ## lb <= x <= ub
4d917a6a858b doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
57 ## A_lb <= A_in*x <= A_ub
9051
1bf0ce0930be Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents: 8920
diff changeset
58 ## @end group
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
59 ## @end example
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10793
diff changeset
60 ##
6741
00116015904d [project @ 2007-06-18 16:07:14 by jwe]
jwe
parents: 6526
diff changeset
61 ## @end ifnottex
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
62 ## @noindent
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
63 ## using a null-space active-set method.
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
64 ##
20165
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
65 ## Any bound (@var{A}, @var{b}, @var{lb}, @var{ub}, @var{A_lb}, @var{A_ub})
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
66 ## may be set to the empty matrix (@code{[]}) if not present. If the initial
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
67 ## guess is feasible the algorithm is faster.
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
68 ##
10059
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
69 ## @table @var
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
70 ## @item options
20165
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
71 ## An optional structure containing the following parameter(s) used to define
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
72 ## the behavior of the solver. Missing elements in the structure take on
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
73 ## default values, so you only need to set the elements that you wish to
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
74 ## change from the default.
10059
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
75 ##
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
76 ## @table @code
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
77 ## @item MaxIter (default: 200)
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
78 ## Maximum number of iterations.
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
79 ## @end table
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
80 ## @end table
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
81 ##
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10793
diff changeset
82 ## @table @var
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10793
diff changeset
83 ## @item info
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10793
diff changeset
84 ## Structure containing run-time information about the algorithm. The
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10793
diff changeset
85 ## following fields are defined:
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10793
diff changeset
86 ##
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
87 ## @table @code
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
88 ## @item solveiter
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
89 ## The number of iterations required to find the solution.
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10793
diff changeset
90 ##
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
91 ## @item info
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10793
diff changeset
92 ## An integer indicating the status of the solution.
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
93 ##
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
94 ## @table @asis
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
95 ## @item 0
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
96 ## The problem is feasible and convex. Global solution found.
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10793
diff changeset
97 ##
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
98 ## @item 1
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
99 ## The problem is not convex. Local solution found.
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10793
diff changeset
100 ##
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
101 ## @item 2
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
102 ## The problem is not convex and unbounded.
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10793
diff changeset
103 ##
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
104 ## @item 3
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
105 ## Maximum number of iterations reached.
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10793
diff changeset
106 ##
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
107 ## @item 6
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
108 ## The problem is infeasible.
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
109 ## @end table
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
110 ## @end table
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10793
diff changeset
111 ## @end table
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
112 ## @end deftypefn
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
113
13027
b9a89ca0fb75 prevent optimization functions from setting ans in workspace at startup
John W. Eaton <jwe@octave.org>
parents: 11589
diff changeset
114 ## PKG_ADD: ## Discard result to avoid polluting workspace with ans at startup.
b9a89ca0fb75 prevent optimization functions from setting ans in workspace at startup
John W. Eaton <jwe@octave.org>
parents: 11589
diff changeset
115 ## PKG_ADD: [~] = __all_opts__ ("qp");
10060
8f51a90eb8d1 implement default opts query and register opts for qp
Jaroslav Hajek <highegg@gmail.com>
parents: 10059
diff changeset
116
10059
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
117 function [x, obj, INFO, lambda] = qp (x0, H, varargin)
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
118
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
119 nargs = nargin;
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
120
10060
8f51a90eb8d1 implement default opts query and register opts for qp
Jaroslav Hajek <highegg@gmail.com>
parents: 10059
diff changeset
121 if (nargin == 1 && ischar (x0) && strcmp (x0, 'defaults'))
8f51a90eb8d1 implement default opts query and register opts for qp
Jaroslav Hajek <highegg@gmail.com>
parents: 10059
diff changeset
122 x = optimset ("MaxIter", 200);
8f51a90eb8d1 implement default opts query and register opts for qp
Jaroslav Hajek <highegg@gmail.com>
parents: 10059
diff changeset
123 return;
8f51a90eb8d1 implement default opts query and register opts for qp
Jaroslav Hajek <highegg@gmail.com>
parents: 10059
diff changeset
124 endif
8f51a90eb8d1 implement default opts query and register opts for qp
Jaroslav Hajek <highegg@gmail.com>
parents: 10059
diff changeset
125
10059
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
126 if (nargs > 2 && isstruct (varargin{end}))
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
127 options = varargin{end};
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
128 nargs--;
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
129 else
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
130 options = struct ();
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
131 endif
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
132
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
133 if (nargs >= 3)
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
134 q = varargin{1};
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
135 else
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
136 q = [];
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
137 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
138
10059
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
139 if (nargs >= 5)
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
140 A = varargin{2};
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
141 b = varargin{3};
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
142 else
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
143 A = [];
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
144 b = [];
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
145 endif
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
146
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
147 if (nargs >= 7)
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
148 lb = varargin{4};
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
149 ub = varargin{5};
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
150 else
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
151 lb = [];
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
152 ub = [];
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
153 endif
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
154
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
155 if (nargs == 10)
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
156 A_lb = varargin{6};
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
157 A_in = varargin{7};
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
158 A_ub = varargin{8};
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
159 else
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
160 A_lb = [];
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
161 A_in = [];
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
162 A_ub = [];
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
163 endif
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
164
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
165 if (nargs == 2 || nargs == 3 || nargs == 5 || nargs == 7 || nargs == 10)
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
166
10064
17ce2a700a97 qp.m: Add missing semicolon.
Ben Abbott <bpabbott@mac.com>
parents: 10060
diff changeset
167 maxit = optimget (options, "MaxIter", 200);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
168
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
169 ## Checking the quadratic penalty
9872
72d6e0de76c7 fix qp, condest and krylov
Jaroslav Hajek <highegg@gmail.com>
parents: 9211
diff changeset
170 if (! issquare (H))
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
171 error ("qp: quadratic penalty matrix not square");
9872
72d6e0de76c7 fix qp, condest and krylov
Jaroslav Hajek <highegg@gmail.com>
parents: 9211
diff changeset
172 elseif (! ishermitian (H))
72d6e0de76c7 fix qp, condest and krylov
Jaroslav Hajek <highegg@gmail.com>
parents: 9211
diff changeset
173 ## warning ("qp: quadratic penalty matrix not hermitian");
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
174 H = (H + H')/2;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
175 endif
9872
72d6e0de76c7 fix qp, condest and krylov
Jaroslav Hajek <highegg@gmail.com>
parents: 9211
diff changeset
176 n = rows (H);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
177
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
178 ## Checking the initial guess (if empty it is resized to the
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
179 ## right dimension and filled with 0)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
180 if (isempty (x0))
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
181 x0 = zeros (n, 1);
10059
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
182 elseif (numel (x0) != n)
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
183 error ("qp: the initial guess has incorrect length");
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
184 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
185
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
186 ## Linear penalty.
10059
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
187 if (isempty (q))
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
188 q = zeros (n, 1);
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
189 elseif (numel (q) != n)
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
190 error ("qp: the linear term has incorrect length");
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
191 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
192
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
193 ## Equality constraint matrices
10059
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
194 if (isempty (A) || isempty (b))
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
195 A = zeros (0, n);
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
196 b = zeros (0, 1);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
197 n_eq = 0;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
198 else
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
199 [n_eq, n1] = size (A);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
200 if (n1 != n)
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
201 error ("qp: equality constraint matrix has incorrect column dimension");
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
202 endif
10059
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
203 if (numel (b) != n_eq)
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
204 error ("qp: equality constraint matrix and vector have inconsistent dimension");
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
205 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
206 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
207
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
208 ## Bound constraints
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
209 Ain = zeros (0, n);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
210 bin = zeros (0, 1);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
211 n_in = 0;
10059
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
212 if (nargs > 5)
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
213 if (! isempty (lb))
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
214 if (numel (lb) != n)
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
215 error ("qp: lower bound has incorrect length");
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
216 elseif (isempty (ub))
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
217 Ain = [Ain; eye(n)];
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
218 bin = [bin; lb];
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
219 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
220 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
221
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
222 if (! isempty (ub))
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
223 if (numel (ub) != n)
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
224 error ("qp: upper bound has incorrect length");
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
225 elseif (isempty (lb))
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
226 Ain = [Ain; -eye(n)];
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
227 bin = [bin; -ub];
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
228 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
229 endif
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
230
8280
5ee11a81688e qp.m: convert b <= x <= b and b <= A*x <= b to equality constraints
Gabriele Pannocchia <g.pannocchia@ing.unipi.it>
parents: 7795
diff changeset
231 if (! isempty (lb) && ! isempty (ub))
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
232 rtol = sqrt (eps);
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
233 for i = 1:n
14868
5d3a684236b0 maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents: 14327
diff changeset
234 if (abs (lb (i) - ub(i)) < rtol*(1 + max (abs (lb(i) + ub(i)))))
8280
5ee11a81688e qp.m: convert b <= x <= b and b <= A*x <= b to equality constraints
Gabriele Pannocchia <g.pannocchia@ing.unipi.it>
parents: 7795
diff changeset
235 ## These are actually an equality constraint
14868
5d3a684236b0 maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents: 14327
diff changeset
236 tmprow = zeros (1,n);
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
237 tmprow(i) = 1;
8280
5ee11a81688e qp.m: convert b <= x <= b and b <= A*x <= b to equality constraints
Gabriele Pannocchia <g.pannocchia@ing.unipi.it>
parents: 7795
diff changeset
238 A = [A;tmprow];
5ee11a81688e qp.m: convert b <= x <= b and b <= A*x <= b to equality constraints
Gabriele Pannocchia <g.pannocchia@ing.unipi.it>
parents: 7795
diff changeset
239 b = [b; 0.5*(lb(i) + ub(i))];
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
240 n_eq = n_eq + 1;
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
241 else
14868
5d3a684236b0 maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents: 14327
diff changeset
242 tmprow = zeros (1,n);
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
243 tmprow(i) = 1;
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
244 Ain = [Ain; tmprow; -tmprow];
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
245 bin = [bin; lb(i); -ub(i)];
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
246 n_in = n_in + 2;
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
247 endif
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
248 endfor
8280
5ee11a81688e qp.m: convert b <= x <= b and b <= A*x <= b to equality constraints
Gabriele Pannocchia <g.pannocchia@ing.unipi.it>
parents: 7795
diff changeset
249 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
250 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
251
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
252 ## Inequality constraints
10059
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
253 if (nargs > 7)
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
254 [dimA_in, n1] = size (A_in);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
255 if (n1 != n)
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
256 error ("qp: inequality constraint matrix has incorrect column dimension");
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
257 else
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
258 if (! isempty (A_lb))
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
259 if (numel (A_lb) != dimA_in)
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
260 error ("qp: inequality constraint matrix and lower bound vector inconsistent");
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
261 elseif (isempty (A_ub))
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
262 Ain = [Ain; A_in];
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
263 bin = [bin; A_lb];
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
264 endif
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
265 endif
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
266 if (! isempty (A_ub))
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
267 if (numel (A_ub) != dimA_in)
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
268 error ("qp: inequality constraint matrix and upper bound vector inconsistent");
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
269 elseif (isempty (A_lb))
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
270 Ain = [Ain; -A_in];
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
271 bin = [bin; -A_ub];
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
272 endif
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
273 endif
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
274
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
275 if (! isempty (A_lb) && ! isempty (A_ub))
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
276 rtol = sqrt (eps);
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
277 for i = 1:dimA_in
19833
9fc020886ae9 maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19697
diff changeset
278 if (abs (A_lb(i) - A_ub(i))
9fc020886ae9 maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19697
diff changeset
279 < rtol*(1 + max (abs (A_lb(i) + A_ub(i)))))
8280
5ee11a81688e qp.m: convert b <= x <= b and b <= A*x <= b to equality constraints
Gabriele Pannocchia <g.pannocchia@ing.unipi.it>
parents: 7795
diff changeset
280 ## These are actually an equality constraint
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
281 tmprow = A_in(i,:);
8280
5ee11a81688e qp.m: convert b <= x <= b and b <= A*x <= b to equality constraints
Gabriele Pannocchia <g.pannocchia@ing.unipi.it>
parents: 7795
diff changeset
282 A = [A;tmprow];
5ee11a81688e qp.m: convert b <= x <= b and b <= A*x <= b to equality constraints
Gabriele Pannocchia <g.pannocchia@ing.unipi.it>
parents: 7795
diff changeset
283 b = [b; 0.5*(A_lb(i) + A_ub(i))];
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
284 n_eq = n_eq + 1;
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
285 else
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
286 tmprow = A_in(i,:);
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
287 Ain = [Ain; tmprow; -tmprow];
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
288 bin = [bin; A_lb(i); -A_ub(i)];
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
289 n_in = n_in + 2;
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
290 endif
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
291 endfor
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
292 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
293 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
294 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
295
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
296 ## Now we should have the following QP:
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
297 ##
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
298 ## min_x 0.5*x'*H*x + x'*q
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
299 ## s.t. A*x = b
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
300 ## Ain*x >= bin
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
301
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
302 ## Discard inequality constraints that have -Inf bounds since those
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
303 ## will never be active.
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
304 idx = isinf (bin) & bin < 0;
6523
589c8dbba916 [project @ 2007-04-12 19:06:44 by jwe]
jwe
parents: 6417
diff changeset
305
6526
2a715c6409a5 [project @ 2007-04-13 20:31:50 by jwe]
jwe
parents: 6523
diff changeset
306 bin(idx) = [];
2a715c6409a5 [project @ 2007-04-13 20:31:50 by jwe]
jwe
parents: 6523
diff changeset
307 Ain(idx,:) = [];
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
308
10059
665ad34efeed qp.m: handle optimset options
Joshua Redstone <redstone@gmail.com>, John W. Eaton <jwe@octave.org>
parents: 9872
diff changeset
309 n_in = numel (bin);
5310
2fbcdc356fc7 [project @ 2005-04-27 15:38:53 by jwe]
jwe
parents: 5307
diff changeset
310
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
311 ## Check if the initial guess is feasible.
19833
9fc020886ae9 maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19697
diff changeset
312 if (isa (x0, "single") || isa (H, "single") || isa (q, "single")
9fc020886ae9 maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19697
diff changeset
313 || isa (A, "single") || isa (b, "single"))
7795
df9519e9990c Handle single precision eps values
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
314 rtol = sqrt (eps ("single"));
df9519e9990c Handle single precision eps values
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
315 else
df9519e9990c Handle single precision eps values
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
316 rtol = sqrt (eps);
df9519e9990c Handle single precision eps values
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
317 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
318
8280
5ee11a81688e qp.m: convert b <= x <= b and b <= A*x <= b to equality constraints
Gabriele Pannocchia <g.pannocchia@ing.unipi.it>
parents: 7795
diff changeset
319 eq_infeasible = (n_eq > 0 && norm (A*x0-b) > rtol*(1+abs (b)));
5ee11a81688e qp.m: convert b <= x <= b and b <= A*x <= b to equality constraints
Gabriele Pannocchia <g.pannocchia@ing.unipi.it>
parents: 7795
diff changeset
320 in_infeasible = (n_in > 0 && any (Ain*x0-bin < -rtol*(1+abs (bin))));
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
321
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
322 info = 0;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
323 if (eq_infeasible || in_infeasible)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
324 ## The initial guess is not feasible.
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
325 ## First define xbar that is feasible with respect to the equality
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
326 ## constraints.
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
327 if (eq_infeasible)
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
328 if (rank (A) < n_eq)
11589
b0084095098e missing semicolons in script files
John W. Eaton <jwe@octave.org>
parents: 11587
diff changeset
329 error ("qp: equality constraint matrix must be full row rank");
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
330 endif
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
331 xbar = pinv (A) * b;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
332 else
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
333 xbar = x0;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
334 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
335
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
336 ## Check if xbar is feasible with respect to the inequality
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
337 ## constraints also.
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
338 if (n_in > 0)
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
339 res = Ain * xbar - bin;
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
340 if (any (res < -rtol * (1 + abs (bin))))
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
341 ## xbar is not feasible with respect to the inequality
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
342 ## constraints. Compute a step in the null space of the
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
343 ## equality constraints, by solving a QP. If the slack is
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
344 ## small, we have a feasible initial guess. Otherwise, the
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
345 ## problem is infeasible.
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
346 if (n_eq > 0)
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
347 Z = null (A);
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
348 if (isempty (Z))
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
349 ## The problem is infeasible because A is square and full
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
350 ## rank, but xbar is not feasible.
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
351 info = 6;
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
352 endif
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
353 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
354
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
355 if (info != 6)
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
356 ## Solve an LP with additional slack variables to find
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
357 ## a feasible starting point.
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
358 gamma = eye (n_in);
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
359 if (n_eq > 0)
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
360 Atmp = [Ain*Z, gamma];
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
361 btmp = -res;
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
362 else
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
363 Atmp = [Ain, gamma];
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
364 btmp = bin;
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
365 endif
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
366 ctmp = [zeros(n-n_eq, 1); ones(n_in, 1)];
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
367 lb = [-Inf(n-n_eq,1); zeros(n_in,1)];
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
368 ub = [];
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
369 ctype = repmat ("L", n_in, 1);
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
370 [P, dummy, status] = glpk (ctmp, Atmp, btmp, lb, ub, ctype);
17897
185038fe7a16 qp.m: Fix test on GLPK's return status (bug #40536)
Julien Bect <julien.bect@supelec.fr>
parents: 17744
diff changeset
371 if ((status == 0)
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
372 && all (abs (P(n-n_eq+1:end)) < rtol * (1 + norm (btmp))))
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
373 ## We found a feasible starting point
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
374 if (n_eq > 0)
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
375 x0 = xbar + Z*P(1:n-n_eq);
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
376 else
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
377 x0 = P(1:n);
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
378 endif
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
379 else
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
380 ## The problem is infeasible
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
381 info = 6;
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
382 endif
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
383 endif
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
384 else
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
385 ## xbar is feasible. We use it a starting point.
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
386 x0 = xbar;
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
387 endif
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
388 else
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
389 ## xbar is feasible. We use it a starting point.
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10540
diff changeset
390 x0 = xbar;
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
391 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
392 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
393
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
394 if (info == 0)
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
395 ## The initial (or computed) guess is feasible.
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
396 ## We call the solver.
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
397 [x, lambda, info, iter] = __qp__ (x0, H, q, A, b, Ain, bin, maxit);
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
398 else
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
399 iter = 0;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
400 x = x0;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
401 lambda = [];
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
402 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
403 obj = 0.5 * x' * H * x + q' * x;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
404 INFO.solveiter = iter;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
405 INFO.info = info;
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
406
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
407 else
6046
34f96dd5441b [project @ 2006-10-10 16:10:25 by jwe]
jwe
parents: 5310
diff changeset
408 print_usage ();
5289
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
409 endif
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
410
b98bf1d70a0a [project @ 2005-04-21 14:43:14 by jwe]
jwe
parents:
diff changeset
411 endfunction
17338
1c89599167a6 maint: End m-files with 1 blank line.
Rik <rik@octave.org>
parents: 14868
diff changeset
412
17897
185038fe7a16 qp.m: Fix test on GLPK's return status (bug #40536)
Julien Bect <julien.bect@supelec.fr>
parents: 17744
diff changeset
413
19833
9fc020886ae9 maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19697
diff changeset
414 ## Test infeasible initial guess (bug #40536)
17931
a404853d2073 qp.m: Condition test on HAVE_GLPK
Mike Miller <mtmiller@ieee.org>
parents: 17897
diff changeset
415 %!testif HAVE_GLPK
17897
185038fe7a16 qp.m: Fix test on GLPK's return status (bug #40536)
Julien Bect <julien.bect@supelec.fr>
parents: 17744
diff changeset
416 %!
185038fe7a16 qp.m: Fix test on GLPK's return status (bug #40536)
Julien Bect <julien.bect@supelec.fr>
parents: 17744
diff changeset
417 %! H = 1; q = 0; # objective: x -> 0.5 x^2
185038fe7a16 qp.m: Fix test on GLPK's return status (bug #40536)
Julien Bect <julien.bect@supelec.fr>
parents: 17744
diff changeset
418 %! A = 1; lb = 1; ub = +inf; # constraint: x >= 1
185038fe7a16 qp.m: Fix test on GLPK's return status (bug #40536)
Julien Bect <julien.bect@supelec.fr>
parents: 17744
diff changeset
419 %! x0 = 0; # infeasible initial guess
185038fe7a16 qp.m: Fix test on GLPK's return status (bug #40536)
Julien Bect <julien.bect@supelec.fr>
parents: 17744
diff changeset
420 %!
185038fe7a16 qp.m: Fix test on GLPK's return status (bug #40536)
Julien Bect <julien.bect@supelec.fr>
parents: 17744
diff changeset
421 %! [x, obj_qp, INFO, lambda] = qp (x0, H, q, [], [], [], [], lb, A, ub);
185038fe7a16 qp.m: Fix test on GLPK's return status (bug #40536)
Julien Bect <julien.bect@supelec.fr>
parents: 17744
diff changeset
422 %!
185038fe7a16 qp.m: Fix test on GLPK's return status (bug #40536)
Julien Bect <julien.bect@supelec.fr>
parents: 17744
diff changeset
423 %! assert (isstruct (INFO) && isfield (INFO, "info") && (INFO.info == 0));
185038fe7a16 qp.m: Fix test on GLPK's return status (bug #40536)
Julien Bect <julien.bect@supelec.fr>
parents: 17744
diff changeset
424 %! assert ([x obj_qp], [1.0 0.5], eps);
185038fe7a16 qp.m: Fix test on GLPK's return status (bug #40536)
Julien Bect <julien.bect@supelec.fr>
parents: 17744
diff changeset
425