File: lsqnonneg.m

package info (click to toggle)
octave 3.8.2-4
  • links: PTS, VCS
  • area: main
  • in suites: jessie, jessie-kfreebsd
  • size: 84,396 kB
  • ctags: 45,547
  • sloc: cpp: 293,356; ansic: 42,041; fortran: 23,669; sh: 13,629; objc: 7,890; yacc: 7,093; lex: 3,442; java: 2,125; makefile: 1,589; perl: 1,009; awk: 974; xml: 34
file content (219 lines) | stat: -rw-r--r-- 6,560 bytes parent folder | download | duplicates (3)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
## Copyright (C) 2008-2013 Bill Denney
## Copyright (C) 2008 Jaroslav Hajek
## Copyright (C) 2009 VZLU Prague
##
## This file is part of Octave.
##
## Octave is free software; you can redistribute it and/or modify it
## under the terms of the GNU General Public License as published by
## the Free Software Foundation; either version 3 of the License, or (at
## your option) any later version.
##
## Octave is distributed in the hope that it will be useful, but
## WITHOUT ANY WARRANTY; without even the implied warranty of
## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
## General Public License for more details.
##
## You should have received a copy of the GNU General Public License
## along with Octave; see the file COPYING.  If not, see
## <http://www.gnu.org/licenses/>.

## -*- texinfo -*-
## @deftypefn  {Function File} {@var{x} =} lsqnonneg (@var{c}, @var{d})
## @deftypefnx {Function File} {@var{x} =} lsqnonneg (@var{c}, @var{d}, @var{x0})
## @deftypefnx {Function File} {@var{x} =} lsqnonneg (@var{c}, @var{d}, @var{x0}, @var{options})
## @deftypefnx {Function File} {[@var{x}, @var{resnorm}] =} lsqnonneg (@dots{})
## @deftypefnx {Function File} {[@var{x}, @var{resnorm}, @var{residual}] =} lsqnonneg (@dots{})
## @deftypefnx {Function File} {[@var{x}, @var{resnorm}, @var{residual}, @var{exitflag}] =} lsqnonneg (@dots{})
## @deftypefnx {Function File} {[@var{x}, @var{resnorm}, @var{residual}, @var{exitflag}, @var{output}] =} lsqnonneg (@dots{})
## @deftypefnx {Function File} {[@var{x}, @var{resnorm}, @var{residual}, @var{exitflag}, @var{output}, @var{lambda}] =} lsqnonneg (@dots{})
## Minimize @code{norm (@var{c}*@var{x} - d)} subject to
## @code{@var{x} >= 0}.  @var{c} and @var{d} must be real.  @var{x0} is an
## optional initial guess for @var{x}.
## Currently, @code{lsqnonneg}
## recognizes these options: @qcode{"MaxIter"}, @qcode{"TolX"}.
## For a description of these options, see @ref{XREFoptimset,,optimset}.
##
## Outputs:
##
## @itemize @bullet
## @item resnorm
##
## The squared 2-norm of the residual: norm (@var{c}*@var{x}-@var{d})^2
##
## @item residual
##
## The residual: @var{d}-@var{c}*@var{x}
##
## @item exitflag
##
## An indicator of convergence.  0 indicates that the iteration count
## was exceeded, and therefore convergence was not reached; >0 indicates
## that the algorithm converged.  (The algorithm is stable and will
## converge given enough iterations.)
##
## @item output
##
## A structure with two fields:
##
## @itemize @bullet
## @item @qcode{"algorithm"}: The algorithm used (@qcode{"nnls"})
##
## @item @qcode{"iterations"}: The number of iterations taken.
## @end itemize
##
## @item lambda
##
## Not implemented.
## @end itemize
## @seealso{optimset, pqpnonneg}
## @end deftypefn

## PKG_ADD: ## Discard result to avoid polluting workspace with ans at startup.
## PKG_ADD: [~] = __all_opts__ ("lsqnonneg");

## This is implemented from Lawson and Hanson's 1973 algorithm on page
## 161 of Solving Least Squares Problems.

function [x, resnorm, residual, exitflag, output, lambda] = lsqnonneg (c, d, x = [], options = struct ())

  if (nargin == 1 && ischar (c) && strcmp (c, 'defaults'))
    x = optimset ("MaxIter", 1e5);
    return;
  endif

  if (! (nargin >= 2 && nargin <= 4 && ismatrix (c) && ismatrix (d) && isstruct (options)))
    print_usage ();
  endif

  ## Lawson-Hanson Step 1 (LH1): initialize the variables.
  m = rows (c);
  n = columns (c);
  if (isempty (x))
    ## Initial guess is 0s.
    x = zeros (n, 1);
  else
    ## ensure nonnegative guess.
    x = max (x, 0);
  endif

  useqr = m >= n;
  max_iter = optimget (options, "MaxIter", 1e5);

  ## Initialize P, according to zero pattern of x.
  p = find (x > 0).';
  if (useqr)
    ## Initialize the QR factorization, economized form.
    [q, r] = qr (c(:,p), 0);
  endif

  iter = 0;

  ## LH3: test for completion.
  while (iter < max_iter)
    while (iter < max_iter)
      iter++;

      ## LH6: compute the positive matrix and find the min norm solution
      ## of the positive problem.
      if (useqr)
        xtmp = r \ q'*d;
      else
        xtmp = c(:,p) \ d;
      endif
      idx = find (xtmp < 0);

      if (isempty (idx))
        ## LH7: tmp solution found, iterate.
        x(:) = 0;
        x(p) = xtmp;
        break;
      else
        ## LH8, LH9: find the scaling factor.
        pidx = p(idx);
        sf = x(pidx)./(x(pidx) - xtmp(idx));
        alpha = min (sf);
        ## LH10: adjust X.
        xx = zeros (n, 1);
        xx(p) = xtmp;
        x += alpha*(xx - x);
        ## LH11: move from P to Z all X == 0.
        ## This corresponds to those indices where minimum of sf is attained.
        idx = idx (sf == alpha);
        p(idx) = [];
        if (useqr)
          ## update the QR factorization.
          [q, r] = qrdelete (q, r, idx);
        endif
      endif
    endwhile

    ## compute the gradient.
    w = c'*(d - c*x);
    w(p) = [];
    tolx = optimget (options, "TolX", 10*eps*norm (c, 1)*length (c));
    if (! any (w > tolx))
      if (useqr)
        ## verify the solution achieved using qr updating.
        ## in the best case, this should only take a single step.
        useqr = false;
        continue;
      else
        ## we're finished.
        break;
      endif
    endif

    ## find the maximum gradient.
    idx = find (w == max (w));
    if (numel (idx) > 1)
      warning ("lsqnonneg:nonunique",
               "a non-unique solution may be returned due to equal gradients");
      idx = idx(1);
    endif
    ## move the index from Z to P. Keep P sorted.
    z = [1:n]; z(p) = [];
    zidx = z(idx);
    jdx = 1 + lookup (p, zidx);
    p = [p(1:jdx-1), zidx, p(jdx:end)];
    if (useqr)
      ## insert the column into the QR factorization.
      [q, r] = qrinsert (q, r, jdx, c(:,zidx));
    endif

  endwhile
  ## LH12: complete.

  ## Generate the additional output arguments.
  if (nargout > 1)
    resnorm = norm (c*x - d) ^ 2;
  endif
  if (nargout > 2)
    residual = d - c*x;
  endif
  exitflag = iter;
  if (nargout > 3 && iter >= max_iter)
    exitflag = 0;
  endif
  if (nargout > 4)
    output = struct ("algorithm", "nnls", "iterations", iter);
  endif
  if (nargout > 5)
    lambda = zeros (size (x));
    lambda(p) = w;
  endif

endfunction


%!test
%! C = [1 0;0 1;2 1];
%! d = [1;3;-2];
%! assert (lsqnonneg (C, d), [0;0.5], 100*eps);

%!test
%! C = [0.0372 0.2869;0.6861 0.7071;0.6233 0.6245;0.6344 0.6170];
%! d = [0.8587;0.1781;0.0747;0.8405];
%! xnew = [0;0.6929];
%! assert (lsqnonneg (C, d), xnew, 0.0001);