GNU Octave  4.0.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
qr.cc
Go to the documentation of this file.
1 /*
2 
3 Copyright (C) 1996-2015 John W. Eaton
4 Copyright (C) 2008-2009 Jaroslav Hajek
5 Copyright (C) 2008-2009 VZLU Prague
6 
7 This file is part of Octave.
8 
9 Octave is free software; you can redistribute it and/or modify it
10 under the terms of the GNU General Public License as published by the
11 Free Software Foundation; either version 3 of the License, or (at your
12 option) any later version.
13 
14 Octave is distributed in the hope that it will be useful, but WITHOUT
15 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
16 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
17 for more details.
18 
19 You should have received a copy of the GNU General Public License
20 along with Octave; see the file COPYING. If not, see
21 <http://www.gnu.org/licenses/>.
22 
23 */
24 
25 #ifdef HAVE_CONFIG_H
26 #include <config.h>
27 #endif
28 
29 #include "CmplxQR.h"
30 #include "CmplxQRP.h"
31 #include "dbleQR.h"
32 #include "dbleQRP.h"
33 #include "fCmplxQR.h"
34 #include "fCmplxQRP.h"
35 #include "floatQR.h"
36 #include "floatQRP.h"
37 #include "SparseQR.h"
38 #include "SparseCmplxQR.h"
39 
40 
41 #include "defun-dld.h"
42 #include "error.h"
43 #include "gripes.h"
44 #include "oct-obj.h"
45 #include "utils.h"
46 
47 template <class MT>
48 static octave_value
49 get_qr_r (const base_qr<MT>& fact)
50 {
51  MT R = fact.R ();
52  if (R.is_square () && fact.regular ())
54  else
55  return R;
56 }
57 
58 // [Q, R] = qr (X): form Q unitary and R upper triangular such
59 // that Q * R = X
60 //
61 // [Q, R] = qr (X, 0): form the economy decomposition such that if X is
62 // m by n then only the first n columns of Q are
63 // computed.
64 //
65 // [Q, R, P] = qr (X): form QRP factorization of X where
66 // P is a permutation matrix such that
67 // A * P = Q * R
68 //
69 // [Q, R, P] = qr (X, 0): form the economy decomposition with
70 // permutation vector P such that Q * R = X (:, P)
71 //
72 // qr (X) alone returns the output of the LAPACK routine dgeqrf, such
73 // that R = triu (qr (X))
74 
75 DEFUN_DLD (qr, args, nargout,
76  "-*- texinfo -*-\n\
77 @deftypefn {Loadable Function} {[@var{Q}, @var{R}, @var{P}] =} qr (@var{A})\n\
78 @deftypefnx {Loadable Function} {[@var{Q}, @var{R}, @var{P}] =} qr (@var{A}, '0')\n\
79 @deftypefnx {Loadable Function} {[@var{C}, @var{R}] =} qr (@var{A}, @var{B})\n\
80 @deftypefnx {Loadable Function} {[@var{C}, @var{R}] =} qr (@var{A}, @var{B}, '0')\n\
81 @cindex QR factorization\n\
82 Compute the QR@tie{}factorization of @var{A}, using standard @sc{lapack}\n\
83 subroutines.\n\
84 \n\
85 For example, given the matrix @code{@var{A} = [1, 2; 3, 4]},\n\
86 \n\
87 @example\n\
88 [@var{Q}, @var{R}] = qr (@var{A})\n\
89 @end example\n\
90 \n\
91 @noindent\n\
92 returns\n\
93 \n\
94 @example\n\
95 @group\n\
96 @var{Q} =\n\
97 \n\
98  -0.31623 -0.94868\n\
99  -0.94868 0.31623\n\
100 \n\
101 @var{R} =\n\
102 \n\
103  -3.16228 -4.42719\n\
104  0.00000 -0.63246\n\
105 @end group\n\
106 @end example\n\
107 \n\
108 The @code{qr} factorization has applications in the solution of least\n\
109 squares problems\n\
110 @tex\n\
111 $$\n\
112 \\min_x \\left\\Vert A x - b \\right\\Vert_2\n\
113 $$\n\
114 @end tex\n\
115 @ifnottex\n\
116 \n\
117 @example\n\
118 min norm(A x - b)\n\
119 @end example\n\
120 \n\
121 @end ifnottex\n\
122 for overdetermined systems of equations (i.e.,\n\
123 @tex\n\
124 $A$\n\
125 @end tex\n\
126 @ifnottex\n\
127 @var{A}\n\
128 @end ifnottex\n\
129 is a tall, thin matrix). The QR@tie{}factorization is\n\
130 @tex\n\
131 $QR = A$ where $Q$ is an orthogonal matrix and $R$ is upper triangular.\n\
132 @end tex\n\
133 @ifnottex\n\
134 @code{@var{Q} * @var{R} = @var{A}} where @var{Q} is an orthogonal matrix and\n\
135 @var{R} is upper triangular.\n\
136 @end ifnottex\n\
137 \n\
138 If given a second argument of @qcode{'0'}, @code{qr} returns an economy-sized\n\
139 QR@tie{}factorization, omitting zero rows of @var{R} and the corresponding\n\
140 columns of @var{Q}.\n\
141 \n\
142 If the matrix @var{A} is full, the permuted QR@tie{}factorization\n\
143 @code{[@var{Q}, @var{R}, @var{P}] = qr (@var{A})} forms the\n\
144 QR@tie{}factorization such that the diagonal entries of @var{R} are\n\
145 decreasing in magnitude order. For example, given the matrix\n\
146 @code{a = [1, 2; 3, 4]},\n\
147 \n\
148 @example\n\
149 [@var{Q}, @var{R}, @var{P}] = qr (@var{A})\n\
150 @end example\n\
151 \n\
152 @noindent\n\
153 returns\n\
154 \n\
155 @example\n\
156 @group\n\
157 @var{Q} =\n\
158 \n\
159  -0.44721 -0.89443\n\
160  -0.89443 0.44721\n\
161 \n\
162 @var{R} =\n\
163 \n\
164  -4.47214 -3.13050\n\
165  0.00000 0.44721\n\
166 \n\
167 @var{P} =\n\
168 \n\
169  0 1\n\
170  1 0\n\
171 @end group\n\
172 @end example\n\
173 \n\
174 The permuted @code{qr} factorization\n\
175 @code{[@var{Q}, @var{R}, @var{P}] = qr (@var{A})} factorization allows the\n\
176 construction of an orthogonal basis of @code{span (A)}.\n\
177 \n\
178 If the matrix @var{A} is sparse, then compute the sparse\n\
179 QR@tie{}factorization of @var{A}, using @sc{CSparse}. As the matrix @var{Q}\n\
180 is in general a full matrix, this function returns the @var{Q}-less\n\
181 factorization @var{R} of @var{A}, such that\n\
182 @code{@var{R} = chol (@var{A}' * @var{A})}.\n\
183 \n\
184 If the final argument is the scalar @code{0} and the number of rows is\n\
185 larger than the number of columns, then an economy factorization is\n\
186 returned. That is @var{R} will have only @code{size (@var{A},1)} rows.\n\
187 \n\
188 If an additional matrix @var{B} is supplied, then @code{qr} returns\n\
189 @var{C}, where @code{@var{C} = @var{Q}' * @var{B}}. This allows the\n\
190 least squares approximation of @code{@var{A} \\ @var{B}} to be calculated\n\
191 as\n\
192 \n\
193 @example\n\
194 @group\n\
195 [@var{C}, @var{R}] = qr (@var{A}, @var{B})\n\
196 x = @var{R} \\ @var{C}\n\
197 @end group\n\
198 @end example\n\
199 @seealso{chol, hess, lu, qz, schur, svd, qrupdate, qrinsert, qrdelete, qrshift}\n\
200 @end deftypefn")
201 {
202  octave_value_list retval;
203 
204  int nargin = args.length ();
205 
206  if (nargin < 1 || nargin > (args(0).is_sparse_type () ? 3 : 2))
207  {
208  print_usage ();
209  return retval;
210  }
211 
212  octave_value arg = args(0);
213 
214  int arg_is_empty = empty_arg ("qr", arg.rows (), arg.columns ());
215 
216  if (arg_is_empty < 0)
217  return retval;
218 
219  if (arg.is_sparse_type ())
220  {
221  bool economy = false;
222  bool is_cmplx = false;
223  int have_b = 0;
224 
225  if (arg.is_complex_type ())
226  is_cmplx = true;
227  if (nargin > 1)
228  {
229  have_b = 1;
230  if (args(nargin-1).is_scalar_type ())
231  {
232  int val = args(nargin-1).int_value ();
233  if (val == 0)
234  {
235  economy = true;
236  have_b = (nargin > 2 ? 2 : 0);
237  }
238  }
239  if (have_b > 0 && args(have_b).is_complex_type ())
240  is_cmplx = true;
241  }
242 
243  if (!error_state)
244  {
245  if (have_b && nargout < 2)
246  error ("qr: incorrect number of output arguments");
247  else if (is_cmplx)
248  {
250  if (!error_state)
251  {
252  if (have_b > 0)
253  {
254  retval(1) = q.R (economy);
255  retval(0) = q.C (args(have_b).complex_matrix_value ());
256  if (arg.rows () < arg.columns ())
257  warning ("qr: non minimum norm solution for under-determined problem");
258  }
259  else if (nargout > 1)
260  {
261  retval(1) = q.R (economy);
262  retval(0) = q.Q ();
263  }
264  else
265  retval(0) = q.R (economy);
266  }
267  }
268  else
269  {
270  SparseQR q (arg.sparse_matrix_value ());
271  if (!error_state)
272  {
273  if (have_b > 0)
274  {
275  retval(1) = q.R (economy);
276  retval(0) = q.C (args(have_b).matrix_value ());
277  if (args(0).rows () < args(0).columns ())
278  warning ("qr: non minimum norm solution for under-determined problem");
279  }
280  else if (nargout > 1)
281  {
282  retval(1) = q.R (economy);
283  retval(0) = q.Q ();
284  }
285  else
286  retval(0) = q.R (economy);
287  }
288  }
289  }
290  }
291  else
292  {
293  QR::type type = (nargout == 0 || nargout == 1) ? QR::raw
294  : nargin == 2
295  ? QR::economy : QR::std;
296 
297  if (arg.is_single_type ())
298  {
299  if (arg.is_real_type ())
300  {
301  FloatMatrix m = arg.float_matrix_value ();
302 
303  if (! error_state)
304  {
305  switch (nargout)
306  {
307  case 0:
308  case 1:
309  {
310  FloatQR fact (m, type);
311  retval(0) = fact.R ();
312  }
313  break;
314 
315  case 2:
316  {
317  FloatQR fact (m, type);
318  retval(1) = get_qr_r (fact);
319  retval(0) = fact.Q ();
320  }
321  break;
322 
323  default:
324  {
325  FloatQRP fact (m, type);
326  if (type == QR::economy)
327  retval(2) = fact.Pvec ();
328  else
329  retval(2) = fact.P ();
330  retval(1) = get_qr_r (fact);
331  retval(0) = fact.Q ();
332  }
333  break;
334  }
335  }
336  }
337  else if (arg.is_complex_type ())
338  {
340 
341  if (! error_state)
342  {
343  switch (nargout)
344  {
345  case 0:
346  case 1:
347  {
348  FloatComplexQR fact (m, type);
349  retval(0) = fact.R ();
350  }
351  break;
352 
353  case 2:
354  {
355  FloatComplexQR fact (m, type);
356  retval(1) = get_qr_r (fact);
357  retval(0) = fact.Q ();
358  }
359  break;
360 
361  default:
362  {
363  FloatComplexQRP fact (m, type);
364  if (type == QR::economy)
365  retval(2) = fact.Pvec ();
366  else
367  retval(2) = fact.P ();
368  retval(1) = get_qr_r (fact);
369  retval(0) = fact.Q ();
370  }
371  break;
372  }
373  }
374  }
375  }
376  else
377  {
378  if (arg.is_real_type ())
379  {
380  Matrix m = arg.matrix_value ();
381 
382  if (! error_state)
383  {
384  switch (nargout)
385  {
386  case 0:
387  case 1:
388  {
389  QR fact (m, type);
390  retval(0) = fact.R ();
391  }
392  break;
393 
394  case 2:
395  {
396  QR fact (m, type);
397  retval(1) = get_qr_r (fact);
398  retval(0) = fact.Q ();
399  }
400  break;
401 
402  default:
403  {
404  QRP fact (m, type);
405  if (type == QR::economy)
406  retval(2) = fact.Pvec ();
407  else
408  retval(2) = fact.P ();
409  retval(1) = get_qr_r (fact);
410  retval(0) = fact.Q ();
411  }
412  break;
413  }
414  }
415  }
416  else if (arg.is_complex_type ())
417  {
419 
420  if (! error_state)
421  {
422  switch (nargout)
423  {
424  case 0:
425  case 1:
426  {
427  ComplexQR fact (m, type);
428  retval(0) = fact.R ();
429  }
430  break;
431 
432  case 2:
433  {
434  ComplexQR fact (m, type);
435  retval(1) = get_qr_r (fact);
436  retval(0) = fact.Q ();
437  }
438  break;
439 
440  default:
441  {
442  ComplexQRP fact (m, type);
443  if (type == QR::economy)
444  retval(2) = fact.Pvec ();
445  else
446  retval(2) = fact.P ();
447  retval(1) = get_qr_r (fact);
448  retval(0) = fact.Q ();
449  }
450  break;
451  }
452  }
453  }
454  else
455  gripe_wrong_type_arg ("qr", arg);
456  }
457  }
458 
459  return retval;
460 }
461 
462 /*
463 %!test
464 %! a = [0, 2, 1; 2, 1, 2];
465 %!
466 %! [q, r] = qr (a);
467 %! [qe, re] = qr (a, 0);
468 %!
469 %! assert (q * r, a, sqrt (eps));
470 %! assert (qe * re, a, sqrt (eps));
471 
472 %!test
473 %! a = [0, 2, 1; 2, 1, 2];
474 %!
475 %! [q, r, p] = qr (a); # FIXME: not giving right dimensions.
476 %! [qe, re, pe] = qr (a, 0);
477 %!
478 %! assert (q * r, a * p, sqrt (eps));
479 %! assert (qe * re, a(:, pe), sqrt (eps));
480 
481 %!test
482 %! a = [0, 2; 2, 1; 1, 2];
483 %!
484 %! [q, r] = qr (a);
485 %! [qe, re] = qr (a, 0);
486 %!
487 %! assert (q * r, a, sqrt (eps));
488 %! assert (qe * re, a, sqrt (eps));
489 
490 %!test
491 %! a = [0, 2; 2, 1; 1, 2];
492 %!
493 %! [q, r, p] = qr (a);
494 %! [qe, re, pe] = qr (a, 0);
495 %!
496 %! assert (q * r, a * p, sqrt (eps));
497 %! assert (qe * re, a(:, pe), sqrt (eps));
498 
499 %!error qr ()
500 %!error qr ([1, 2; 3, 4], 0, 2)
501 
502 %!function retval = __testqr (q, r, a, p)
503 %! tol = 100*eps (class (q));
504 %! retval = 0;
505 %! if (nargin == 3)
506 %! n1 = norm (q*r - a);
507 %! n2 = norm (q'*q - eye (columns (q)));
508 %! retval = (n1 < tol && n2 < tol);
509 %! else
510 %! n1 = norm (q'*q - eye (columns (q)));
511 %! retval = (n1 < tol);
512 %! if (isvector (p))
513 %! n2 = norm (q*r - a(:,p));
514 %! retval = (retval && n2 < tol);
515 %! else
516 %! n2 = norm (q*r - a*p);
517 %! retval = (retval && n2 < tol);
518 %! endif
519 %! endif
520 %!endfunction
521 
522 %!test
523 %! t = ones (24, 1);
524 %! j = 1;
525 %!
526 %! if (false) # eliminate big matrix tests
527 %! a = rand (5000, 20);
528 %! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
529 %! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
530 %! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
531 %! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
532 %!
533 %! a = a+1i*eps;
534 %! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
535 %! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
536 %! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
537 %! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
538 %! endif
539 %!
540 %! a = [ ones(1,15); sqrt(eps)*eye(15) ];
541 %! [q,r] = qr (a); t(j++) = __testqr (q, r, a);
542 %! [q,r] = qr (a'); t(j++) = __testqr (q, r, a');
543 %! [q,r,p] = qr (a); t(j++) = __testqr (q, r, a, p);
544 %! [q,r,p] = qr (a'); t(j++) = __testqr (q, r, a', p);
545 %!
546 %! a = a+1i*eps;
547 %! [q,r] = qr (a); t(j++) = __testqr (q, r, a);
548 %! [q,r] = qr (a'); t(j++) = __testqr (q, r, a');
549 %! [q,r,p] = qr (a); t(j++) = __testqr (q, r, a, p);
550 %! [q,r,p] = qr (a'); t(j++) = __testqr (q, r, a', p);
551 %!
552 %! a = [ ones(1,15); sqrt(eps)*eye(15) ];
553 %! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
554 %! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
555 %! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
556 %! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
557 %!
558 %! a = a+1i*eps;
559 %! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
560 %! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
561 %! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
562 %! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
563 %!
564 %! a = [ 611 196 -192 407 -8 -52 -49 29
565 %! 196 899 113 -192 -71 -43 -8 -44
566 %! -192 113 899 196 61 49 8 52
567 %! 407 -192 196 611 8 44 59 -23
568 %! -8 -71 61 8 411 -599 208 208
569 %! -52 -43 49 44 -599 411 208 208
570 %! -49 -8 8 59 208 208 99 -911
571 %! 29 -44 52 -23 208 208 -911 99 ];
572 %! [q,r] = qr (a);
573 %!
574 %! assert (all (t) && norm (q*r - a) < 5000*eps);
575 
576 %!test
577 %! a = single ([0, 2, 1; 2, 1, 2]);
578 %!
579 %! [q, r] = qr (a);
580 %! [qe, re] = qr (a, 0);
581 %!
582 %! assert (q * r, a, sqrt (eps ("single")));
583 %! assert (qe * re, a, sqrt (eps ("single")));
584 
585 %!test
586 %! a = single ([0, 2, 1; 2, 1, 2]);
587 %!
588 %! [q, r, p] = qr (a); # FIXME: not giving right dimensions.
589 %! [qe, re, pe] = qr (a, 0);
590 %!
591 %! assert (q * r, a * p, sqrt (eps ("single")));
592 %! assert (qe * re, a(:, pe), sqrt (eps ("single")));
593 
594 %!test
595 %! a = single ([0, 2; 2, 1; 1, 2]);
596 %!
597 %! [q, r] = qr (a);
598 %! [qe, re] = qr (a, 0);
599 %!
600 %! assert (q * r, a, sqrt (eps ("single")));
601 %! assert (qe * re, a, sqrt (eps ("single")));
602 
603 %!test
604 %! a = single ([0, 2; 2, 1; 1, 2]);
605 %!
606 %! [q, r, p] = qr (a);
607 %! [qe, re, pe] = qr (a, 0);
608 %!
609 %! assert (q * r, a * p, sqrt (eps ("single")));
610 %! assert (qe * re, a(:, pe), sqrt (eps ("single")));
611 
612 %!error qr ()
613 %!error qr ([1, 2; 3, 4], 0, 2)
614 
615 %!test
616 %! t = ones (24, 1);
617 %! j = 1;
618 %!
619 %! if (false) # eliminate big matrix tests
620 %! a = rand (5000,20);
621 %! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
622 %! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
623 %! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
624 %! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
625 %!
626 %! a = a+1i*eps ("single");
627 %! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
628 %! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
629 %! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
630 %! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
631 %! endif
632 %!
633 %! a = [ ones(1,15); sqrt(eps("single"))*eye(15) ];
634 %! [q,r] = qr (a); t(j++) = __testqr (q, r, a);
635 %! [q,r] = qr (a'); t(j++) = __testqr (q, r, a');
636 %! [q,r,p] = qr (a); t(j++) = __testqr (q, r, a, p);
637 %! [q,r,p] = qr (a'); t(j++) = __testqr (q, r, a', p);
638 %!
639 %! a = a+1i*eps ("single");
640 %! [q,r] = qr (a); t(j++) = __testqr (q, r, a);
641 %! [q,r] = qr (a'); t(j++) = __testqr (q, r, a');
642 %! [q,r,p] = qr (a); t(j++) = __testqr (q, r, a, p);
643 %! [q,r,p] = qr (a'); t(j++) = __testqr (q, r, a', p);
644 %!
645 %! a = [ ones(1,15); sqrt(eps("single"))*eye(15) ];
646 %! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
647 %! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
648 %! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
649 %! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
650 %!
651 %! a = a+1i*eps ("single");
652 %! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
653 %! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
654 %! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
655 %! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a',p);
656 %!
657 %! a = [ 611 196 -192 407 -8 -52 -49 29
658 %! 196 899 113 -192 -71 -43 -8 -44
659 %! -192 113 899 196 61 49 8 52
660 %! 407 -192 196 611 8 44 59 -23
661 %! -8 -71 61 8 411 -599 208 208
662 %! -52 -43 49 44 -599 411 208 208
663 %! -49 -8 8 59 208 208 99 -911
664 %! 29 -44 52 -23 208 208 -911 99 ];
665 %! [q,r] = qr (a);
666 %!
667 %! assert (all (t) && norm (q*r-a) < 5000*eps ("single"));
668 
669 ## The deactivated tests below can't be tested till rectangular back-subs is
670 ## implemented for sparse matrices.
671 
672 %!testif HAVE_CXSPARSE
673 %! n = 20; d = 0.2;
674 %! a = sprandn (n,n,d) + speye (n,n);
675 %! r = qr (a);
676 %! assert (r'*r, a'*a, 1e-10)
677 
678 %!testif HAVE_COLAMD
679 %! n = 20; d = 0.2;
680 %! a = sprandn (n,n,d) + speye (n,n);
681 %! q = symamd (a);
682 %! a = a(q,q);
683 %! r = qr (a);
684 %! assert (r'*r, a'*a, 1e-10)
685 
686 %!testif HAVE_CXSPARSE
687 %! n = 20; d = 0.2;
688 %! a = sprandn (n,n,d) + speye (n,n);
689 %! [c,r] = qr (a, ones (n,1));
690 %! assert (r\c, full (a)\ones (n,1), 10e-10)
691 
692 %!testif HAVE_CXSPARSE
693 %! n = 20; d = 0.2;
694 %! a = sprandn (n,n,d) + speye (n,n);
695 %! b = randn (n,2);
696 %! [c,r] = qr (a, b);
697 %! assert (r\c, full (a)\b, 10e-10)
698 
699 %% Test under-determined systems!!
700 %!#testif HAVE_CXSPARSE
701 %! n = 20; d = 0.2;
702 %! a = sprandn (n,n+1,d) + speye (n,n+1);
703 %! b = randn (n,2);
704 %! [c,r] = qr (a, b);
705 %! assert (r\c, full (a)\b, 10e-10)
706 
707 %!testif HAVE_CXSPARSE
708 %! n = 20; d = 0.2;
709 %! a = 1i*sprandn (n,n,d) + speye (n,n);
710 %! r = qr (a);
711 %! assert (r'*r,a'*a,1e-10)
712 
713 %!testif HAVE_COLAMD
714 %! n = 20; d = 0.2;
715 %! a = 1i*sprandn (n,n,d) + speye (n,n);
716 %! q = symamd (a);
717 %! a = a(q,q);
718 %! r = qr (a);
719 %! assert (r'*r, a'*a, 1e-10)
720 
721 %!testif HAVE_CXSPARSE
722 %! n = 20; d = 0.2;
723 %! a = 1i*sprandn (n,n,d) + speye (n,n);
724 %! [c,r] = qr (a, ones (n,1));
725 %! assert (r\c, full (a)\ones (n,1), 10e-10)
726 
727 %!testif HAVE_CXSPARSE
728 %! n = 20; d = 0.2;
729 %! a = 1i*sprandn (n,n,d) + speye (n,n);
730 %! b = randn (n,2);
731 %! [c,r] = qr (a, b);
732 %! assert (r\c, full (a)\b, 10e-10)
733 
734 %% Test under-determined systems!!
735 %!#testif HAVE_CXSPARSE
736 %! n = 20; d = 0.2;
737 %! a = 1i*sprandn (n,n+1,d) + speye (n,n+1);
738 %! b = randn (n,2);
739 %! [c,r] = qr (a, b);
740 %! assert (r\c, full (a)\b, 10e-10)
741 
742 %!error qr (sprandn (10,10,0.2), ones (10,1))
743 */
744 
745 static
746 bool check_qr_dims (const octave_value& q, const octave_value& r,
747  bool allow_ecf = false)
748 {
749  octave_idx_type m = q.rows ();
750  octave_idx_type k = r.rows ();
751  octave_idx_type n = r.columns ();
752  return ((q.ndims () == 2 && r.ndims () == 2 && k == q.columns ())
753  && (m == k || (allow_ecf && k == n && k < m)));
754 }
755 
756 static
757 bool check_index (const octave_value& i, bool vector_allowed = false)
758 {
759  return ((i.is_real_type () || i.is_integer_type ())
760  && (i.is_scalar_type () || vector_allowed));
761 }
762 
763 DEFUN_DLD (qrupdate, args, ,
764  "-*- texinfo -*-\n\
765 @deftypefn {Loadable Function} {[@var{Q1}, @var{R1}] =} qrupdate (@var{Q}, @var{R}, @var{u}, @var{v})\n\
766 Given a QR@tie{}factorization of a real or complex matrix\n\
767 @w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and\n\
768 @var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization of\n\
769 @w{@var{A} + @var{u}*@var{v}'}, where @var{u} and @var{v} are column vectors\n\
770 (rank-1 update) or matrices with equal number of columns\n\
771 (rank-k update). Notice that the latter case is done as a sequence of rank-1\n\
772 updates; thus, for k large enough, it will be both faster and more accurate\n\
773 to recompute the factorization from scratch.\n\
774 \n\
775 The QR@tie{}factorization supplied may be either full (Q is square) or\n\
776 economized (R is square).\n\
777 \n\
778 @seealso{qr, qrinsert, qrdelete, qrshift}\n\
779 @end deftypefn")
780 {
781  octave_idx_type nargin = args.length ();
782  octave_value_list retval;
783 
784  if (nargin != 4)
785  {
786  print_usage ();
787  return retval;
788  }
789 
790  octave_value argq = args(0);
791  octave_value argr = args(1);
792  octave_value argu = args(2);
793  octave_value argv = args(3);
794 
795  if (argq.is_numeric_type () && argr.is_numeric_type ()
796  && argu.is_numeric_type () && argv.is_numeric_type ())
797  {
798  if (check_qr_dims (argq, argr, true))
799  {
800  if (argq.is_real_type ()
801  && argr.is_real_type ()
802  && argu.is_real_type ()
803  && argv.is_real_type ())
804  {
805  // all real case
806  if (argq.is_single_type ()
807  || argr.is_single_type ()
808  || argu.is_single_type ()
809  || argv.is_single_type ())
810  {
811  FloatMatrix Q = argq.float_matrix_value ();
812  FloatMatrix R = argr.float_matrix_value ();
813  FloatMatrix u = argu.float_matrix_value ();
814  FloatMatrix v = argv.float_matrix_value ();
815 
816  FloatQR fact (Q, R);
817  fact.update (u, v);
818 
819  retval(1) = get_qr_r (fact);
820  retval(0) = fact.Q ();
821  }
822  else
823  {
824  Matrix Q = argq.matrix_value ();
825  Matrix R = argr.matrix_value ();
826  Matrix u = argu.matrix_value ();
827  Matrix v = argv.matrix_value ();
828 
829  QR fact (Q, R);
830  fact.update (u, v);
831 
832  retval(1) = get_qr_r (fact);
833  retval(0) = fact.Q ();
834  }
835  }
836  else
837  {
838  // complex case
839  if (argq.is_single_type ()
840  || argr.is_single_type ()
841  || argu.is_single_type ()
842  || argv.is_single_type ())
843  {
848 
849  FloatComplexQR fact (Q, R);
850  fact.update (u, v);
851 
852  retval(1) = get_qr_r (fact);
853  retval(0) = fact.Q ();
854  }
855  else
856  {
861 
862  ComplexQR fact (Q, R);
863  fact.update (u, v);
864 
865  retval(1) = get_qr_r (fact);
866  retval(0) = fact.Q ();
867  }
868  }
869  }
870  else
871  error ("qrupdate: Q and R dimensions don't match");
872  }
873  else
874  error ("qrupdate: Q, R, U, and V must be numeric");
875 
876  return retval;
877 }
878 
879 /*
880 %!shared A, u, v, Ac, uc, vc
881 %! A = [0.091364 0.613038 0.999083;
882 %! 0.594638 0.425302 0.603537;
883 %! 0.383594 0.291238 0.085574;
884 %! 0.265712 0.268003 0.238409;
885 %! 0.669966 0.743851 0.445057 ];
886 %!
887 %! u = [0.85082;
888 %! 0.76426;
889 %! 0.42883;
890 %! 0.53010;
891 %! 0.80683 ];
892 %!
893 %! v = [0.98810;
894 %! 0.24295;
895 %! 0.43167 ];
896 %!
897 %! Ac = [0.620405 + 0.956953i 0.480013 + 0.048806i 0.402627 + 0.338171i;
898 %! 0.589077 + 0.658457i 0.013205 + 0.279323i 0.229284 + 0.721929i;
899 %! 0.092758 + 0.345687i 0.928679 + 0.241052i 0.764536 + 0.832406i;
900 %! 0.912098 + 0.721024i 0.049018 + 0.269452i 0.730029 + 0.796517i;
901 %! 0.112849 + 0.603871i 0.486352 + 0.142337i 0.355646 + 0.151496i ];
902 %!
903 %! uc = [0.20351 + 0.05401i;
904 %! 0.13141 + 0.43708i;
905 %! 0.29808 + 0.08789i;
906 %! 0.69821 + 0.38844i;
907 %! 0.74871 + 0.25821i ];
908 %!
909 %! vc = [0.85839 + 0.29468i;
910 %! 0.20820 + 0.93090i;
911 %! 0.86184 + 0.34689i ];
912 %!
913 
914 %!test
915 %! [Q,R] = qr (A);
916 %! [Q,R] = qrupdate (Q, R, u, v);
917 %! assert (norm (vec (Q'*Q - eye (5)), Inf) < 1e1*eps);
918 %! assert (norm (vec (triu (R)-R), Inf) == 0);
919 %! assert (norm (vec (Q*R - A - u*v'), Inf) < norm (A)*1e1*eps);
920 %!
921 %!test
922 %! [Q,R] = qr (Ac);
923 %! [Q,R] = qrupdate (Q, R, uc, vc);
924 %! assert (norm (vec (Q'*Q - eye (5)), Inf) < 1e1*eps);
925 %! assert (norm (vec (triu (R)-R), Inf) == 0);
926 %! assert (norm (vec (Q*R - Ac - uc*vc'), Inf) < norm (Ac)*1e1*eps);
927 
928 %!test
929 %! [Q,R] = qr (single (A));
930 %! [Q,R] = qrupdate (Q, R, single (u), single (v));
931 %! assert (norm (vec (Q'*Q - eye (5,"single")), Inf) < 1e1*eps ("single"));
932 %! assert (norm (vec (triu (R)-R), Inf) == 0);
933 %! assert (norm (vec (Q*R - single (A) - single (u)*single (v)'), Inf) < norm (single (A))*1e1*eps ("single"));
934 %!
935 %!test
936 %! [Q,R] = qr (single (Ac));
937 %! [Q,R] = qrupdate (Q, R, single (uc), single (vc));
938 %! assert (norm (vec (Q'*Q - eye (5,"single")), Inf) < 1e1*eps ("single"));
939 %! assert (norm (vec (triu (R)-R), Inf) == 0);
940 %! assert (norm (vec (Q*R - single (Ac) - single (uc)*single (vc)'), Inf) < norm (single (Ac))*1e1*eps ("single"));
941 */
942 
943 DEFUN_DLD (qrinsert, args, ,
944  "-*- texinfo -*-\n\
945 @deftypefn {Loadable Function} {[@var{Q1}, @var{R1}] =} qrinsert (@var{Q}, @var{R}, @var{j}, @var{x}, @var{orient})\n\
946 Given a QR@tie{}factorization of a real or complex matrix\n\
947 @w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and\n\
948 @var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization of\n\
949 @w{[A(:,1:j-1) x A(:,j:n)]}, where @var{u} is a column vector to be inserted\n\
950 into @var{A} (if @var{orient} is @qcode{\"col\"}), or the\n\
951 QR@tie{}factorization of @w{[A(1:j-1,:);x;A(:,j:n)]}, where @var{x} is a row\n\
952 vector to be inserted into @var{A} (if @var{orient} is @qcode{\"row\"}).\n\
953 \n\
954 The default value of @var{orient} is @qcode{\"col\"}. If @var{orient} is\n\
955 @qcode{\"col\"}, @var{u} may be a matrix and @var{j} an index vector\n\
956 resulting in the QR@tie{}factorization of a matrix @var{B} such that\n\
957 @w{B(:,@var{j})} gives @var{u} and @w{B(:,@var{j}) = []} gives @var{A}.\n\
958 Notice that the latter case is done as a sequence of k insertions;\n\
959 thus, for k large enough, it will be both faster and more accurate to\n\
960 recompute the factorization from scratch.\n\
961 \n\
962 If @var{orient} is @qcode{\"col\"}, the QR@tie{}factorization supplied may\n\
963 be either full (Q is square) or economized (R is square).\n\
964 \n\
965 If @var{orient} is @qcode{\"row\"}, full factorization is needed.\n\
966 @seealso{qr, qrupdate, qrdelete, qrshift}\n\
967 @end deftypefn")
968 {
969  octave_idx_type nargin = args.length ();
970  octave_value_list retval;
971 
972  if (nargin < 4 || nargin > 5)
973  {
974  print_usage ();
975  return retval;
976  }
977 
978  octave_value argq = args(0);
979  octave_value argr = args(1);
980  octave_value argj = args(2);
981  octave_value argx = args(3);
982 
983  if (argq.is_numeric_type () && argr.is_numeric_type ()
984  && argx.is_numeric_type ()
985  && (nargin < 5 || args(4).is_string ()))
986  {
987  std::string orient = (nargin < 5) ? "col" : args(4).string_value ();
988 
989  bool col = orient == "col";
990 
991  if (col || orient == "row")
992  if (check_qr_dims (argq, argr, col)
993  && (col || argx.rows () == 1))
994  {
995  if (check_index (argj, col))
996  {
999 
1000  octave_idx_type one = 1;
1001 
1002  if (argq.is_real_type ()
1003  && argr.is_real_type ()
1004  && argx.is_real_type ())
1005  {
1006  // real case
1007  if (argq.is_single_type ()
1008  || argr.is_single_type ()
1009  || argx.is_single_type ())
1010  {
1011  FloatMatrix Q = argq.float_matrix_value ();
1012  FloatMatrix R = argr.float_matrix_value ();
1013  FloatMatrix x = argx.float_matrix_value ();
1014 
1015  FloatQR fact (Q, R);
1016 
1017  if (col)
1018  fact.insert_col (x, j-one);
1019  else
1020  fact.insert_row (x.row (0), j(0)-one);
1021 
1022  retval(1) = get_qr_r (fact);
1023  retval(0) = fact.Q ();
1024 
1025  }
1026  else
1027  {
1028  Matrix Q = argq.matrix_value ();
1029  Matrix R = argr.matrix_value ();
1030  Matrix x = argx.matrix_value ();
1031 
1032  QR fact (Q, R);
1033 
1034  if (col)
1035  fact.insert_col (x, j-one);
1036  else
1037  fact.insert_row (x.row (0), j(0)-one);
1038 
1039  retval(1) = get_qr_r (fact);
1040  retval(0) = fact.Q ();
1041 
1042  }
1043  }
1044  else
1045  {
1046  // complex case
1047  if (argq.is_single_type ()
1048  || argr.is_single_type ()
1049  || argx.is_single_type ())
1050  {
1053  FloatComplexMatrix R =
1057 
1058  FloatComplexQR fact (Q, R);
1059 
1060  if (col)
1061  fact.insert_col (x, j-one);
1062  else
1063  fact.insert_row (x.row (0), j(0)-one);
1064 
1065  retval(1) = get_qr_r (fact);
1066  retval(0) = fact.Q ();
1067  }
1068  else
1069  {
1071  ComplexMatrix R = argr.complex_matrix_value ();
1073 
1074  ComplexQR fact (Q, R);
1075 
1076  if (col)
1077  fact.insert_col (x, j-one);
1078  else
1079  fact.insert_row (x.row (0), j(0)-one);
1080 
1081  retval(1) = get_qr_r (fact);
1082  retval(0) = fact.Q ();
1083  }
1084  }
1085 
1086  }
1087  else
1088  error ("qrinsert: invalid index J");
1089  }
1090  else
1091  error ("qrinsert: dimension mismatch");
1092 
1093  else
1094  error ("qrinsert: ORIENT must be \"col\" or \"row\"");
1095  }
1096  else
1097  print_usage ();
1098 
1099  return retval;
1100 }
1101 
1102 /*
1103 %!test
1104 %! [Q,R] = qr (A);
1105 %! [Q,R] = qrinsert (Q, R, 3, u);
1106 %! assert (norm (vec (Q'*Q - eye (5)), Inf) < 1e1*eps);
1107 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1108 %! assert (norm (vec (Q*R - [A(:,1:2) u A(:,3)]), Inf) < norm (A)*1e1*eps);
1109 %!test
1110 %! [Q,R] = qr (Ac);
1111 %! [Q,R] = qrinsert (Q, R, 3, uc);
1112 %! assert (norm (vec (Q'*Q - eye (5)), Inf) < 1e1*eps);
1113 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1114 %! assert (norm (vec (Q*R - [Ac(:,1:2) uc Ac(:,3)]), Inf) < norm (Ac)*1e1*eps);
1115 %!test
1116 %! x = [0.85082 0.76426 0.42883 ];
1117 %!
1118 %! [Q,R] = qr (A);
1119 %! [Q,R] = qrinsert (Q, R, 3, x, "row");
1120 %! assert (norm (vec (Q'*Q - eye (6)), Inf) < 1e1*eps);
1121 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1122 %! assert (norm (vec (Q*R - [A(1:2,:);x;A(3:5,:)]), Inf) < norm (A)*1e1*eps);
1123 %!test
1124 %! x = [0.20351 + 0.05401i 0.13141 + 0.43708i 0.29808 + 0.08789i ];
1125 %!
1126 %! [Q,R] = qr (Ac);
1127 %! [Q,R] = qrinsert (Q, R, 3, x, "row");
1128 %! assert (norm (vec (Q'*Q - eye (6)), Inf) < 1e1*eps);
1129 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1130 %! assert (norm (vec (Q*R - [Ac(1:2,:);x;Ac(3:5,:)]), Inf) < norm (Ac)*1e1*eps);
1131 
1132 %!test
1133 %! [Q,R] = qr (single (A));
1134 %! [Q,R] = qrinsert (Q, R, 3, single (u));
1135 %! assert (norm (vec (Q'*Q - eye (5,"single")), Inf) < 1e1*eps ("single"));
1136 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1137 %! assert (norm (vec (Q*R - single ([A(:,1:2) u A(:,3)])), Inf) < norm (single (A))*1e1*eps ("single"));
1138 %!test
1139 %! [Q,R] = qr (single (Ac));
1140 %! [Q,R] = qrinsert (Q, R, 3, single (uc));
1141 %! assert (norm (vec (Q'*Q - eye (5,"single")), Inf) < 1e1*eps ("single"));
1142 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1143 %! assert (norm (vec (Q*R - single ([Ac(:,1:2) uc Ac(:,3)])), Inf) < norm (single (Ac))*1e1*eps ("single"));
1144 %!test
1145 %! x = single ([0.85082 0.76426 0.42883 ]);
1146 %!
1147 %! [Q,R] = qr (single (A));
1148 %! [Q,R] = qrinsert (Q, R, 3, x, "row");
1149 %! assert (norm (vec (Q'*Q - eye (6,"single")), Inf) < 1e1*eps ("single"));
1150 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1151 %! assert (norm (vec (Q*R - single ([A(1:2,:);x;A(3:5,:)])), Inf) < norm (single (A))*1e1*eps ("single"));
1152 %!test
1153 %! x = single ([0.20351 + 0.05401i 0.13141 + 0.43708i 0.29808 + 0.08789i ]);
1154 %!
1155 %! [Q,R] = qr (single (Ac));
1156 %! [Q,R] = qrinsert (Q, R, 3, x, "row");
1157 %! assert (norm (vec (Q'*Q - eye (6,"single")), Inf) < 1e1*eps ("single"));
1158 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1159 %! assert (norm (vec (Q*R - single ([Ac(1:2,:);x;Ac(3:5,:)])), Inf) < norm (single (Ac))*1e1*eps ("single"));
1160 */
1161 
1162 DEFUN_DLD (qrdelete, args, ,
1163  "-*- texinfo -*-\n\
1164 @deftypefn {Loadable Function} {[@var{Q1}, @var{R1}] =} qrdelete (@var{Q}, @var{R}, @var{j}, @var{orient})\n\
1165 Given a QR@tie{}factorization of a real or complex matrix\n\
1166 @w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and\n\
1167 @var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization of\n\
1168 @w{[A(:,1:j-1) A(:,j+1:n)]}, i.e., @var{A} with one column deleted\n\
1169 (if @var{orient} is @qcode{\"col\"}), or the QR@tie{}factorization of\n\
1170 @w{[A(1:j-1,:);A(j+1:n,:)]}, i.e., @var{A} with one row deleted (if\n\
1171 @var{orient} is @qcode{\"row\"}).\n\
1172 \n\
1173 The default value of @var{orient} is @qcode{\"col\"}.\n\
1174 \n\
1175 If @var{orient} is @qcode{\"col\"}, @var{j} may be an index vector\n\
1176 resulting in the QR@tie{}factorization of a matrix @var{B} such that\n\
1177 @w{A(:,@var{j}) = []} gives @var{B}. Notice that the latter case is done as\n\
1178 a sequence of k deletions; thus, for k large enough, it will be both faster\n\
1179 and more accurate to recompute the factorization from scratch.\n\
1180 \n\
1181 If @var{orient} is @qcode{\"col\"}, the QR@tie{}factorization supplied may\n\
1182 be either full (Q is square) or economized (R is square).\n\
1183 \n\
1184 If @var{orient} is @qcode{\"row\"}, full factorization is needed.\n\
1185 @seealso{qr, qrupdate, qrinsert, qrshift}\n\
1186 @end deftypefn")
1187 {
1188  octave_idx_type nargin = args.length ();
1189  octave_value_list retval;
1190 
1191  if (nargin < 3 || nargin > 4)
1192  {
1193  print_usage ();
1194  return retval;
1195  }
1196 
1197  octave_value argq = args(0);
1198  octave_value argr = args(1);
1199  octave_value argj = args(2);
1200 
1201  if (argq.is_numeric_type () && argr.is_numeric_type ()
1202  && (nargin < 4 || args(3).is_string ()))
1203  {
1204  std::string orient = (nargin < 4) ? "col" : args(3).string_value ();
1205 
1206  bool col = orient == "col";
1207 
1208  if (col || orient == "row")
1209  if (check_qr_dims (argq, argr, col))
1210  {
1211  if (check_index (argj, col))
1212  {
1214  = argj.octave_idx_type_vector_value ();
1215 
1216  octave_idx_type one = 1;
1217 
1218  if (argq.is_real_type ()
1219  && argr.is_real_type ())
1220  {
1221  // real case
1222  if (argq.is_single_type ()
1223  || argr.is_single_type ())
1224  {
1225  FloatMatrix Q = argq.float_matrix_value ();
1226  FloatMatrix R = argr.float_matrix_value ();
1227 
1228  FloatQR fact (Q, R);
1229 
1230  if (col)
1231  fact.delete_col (j-one);
1232  else
1233  fact.delete_row (j(0)-one);
1234 
1235  retval(1) = get_qr_r (fact);
1236  retval(0) = fact.Q ();
1237  }
1238  else
1239  {
1240  Matrix Q = argq.matrix_value ();
1241  Matrix R = argr.matrix_value ();
1242 
1243  QR fact (Q, R);
1244 
1245  if (col)
1246  fact.delete_col (j-one);
1247  else
1248  fact.delete_row (j(0)-one);
1249 
1250  retval(1) = get_qr_r (fact);
1251  retval(0) = fact.Q ();
1252  }
1253  }
1254  else
1255  {
1256  // complex case
1257  if (argq.is_single_type ()
1258  || argr.is_single_type ())
1259  {
1262  FloatComplexMatrix R =
1264 
1265  FloatComplexQR fact (Q, R);
1266 
1267  if (col)
1268  fact.delete_col (j-one);
1269  else
1270  fact.delete_row (j(0)-one);
1271 
1272  retval(1) = get_qr_r (fact);
1273  retval(0) = fact.Q ();
1274  }
1275  else
1276  {
1278  ComplexMatrix R = argr.complex_matrix_value ();
1279 
1280  ComplexQR fact (Q, R);
1281 
1282  if (col)
1283  fact.delete_col (j-one);
1284  else
1285  fact.delete_row (j(0)-one);
1286 
1287  retval(1) = get_qr_r (fact);
1288  retval(0) = fact.Q ();
1289  }
1290  }
1291  }
1292  else
1293  error ("qrdelete: invalid index J");
1294  }
1295  else
1296  error ("qrdelete: dimension mismatch");
1297 
1298  else
1299  error ("qrdelete: ORIENT must be \"col\" or \"row\"");
1300  }
1301  else
1302  print_usage ();
1303 
1304  return retval;
1305 }
1306 
1307 /*
1308 %!test
1309 %! AA = [0.091364 0.613038 0.027504 0.999083;
1310 %! 0.594638 0.425302 0.562834 0.603537;
1311 %! 0.383594 0.291238 0.742073 0.085574;
1312 %! 0.265712 0.268003 0.783553 0.238409;
1313 %! 0.669966 0.743851 0.457255 0.445057 ];
1314 %!
1315 %! [Q,R] = qr (AA);
1316 %! [Q,R] = qrdelete (Q, R, 3);
1317 %! assert (norm (vec (Q'*Q - eye (5)), Inf) < 16*eps);
1318 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1319 %! assert (norm (vec (Q*R - [AA(:,1:2) AA(:,4)]), Inf) < norm (AA)*1e1*eps);
1320 %!
1321 %!test
1322 %! AA = [0.364554 + 0.993117i 0.669818 + 0.510234i 0.426568 + 0.041337i 0.847051 + 0.233291i;
1323 %! 0.049600 + 0.242783i 0.448946 + 0.484022i 0.141155 + 0.074420i 0.446746 + 0.392706i;
1324 %! 0.581922 + 0.657416i 0.581460 + 0.030016i 0.219909 + 0.447288i 0.201144 + 0.069132i;
1325 %! 0.694986 + 0.000571i 0.682327 + 0.841712i 0.807537 + 0.166086i 0.192767 + 0.358098i;
1326 %! 0.945002 + 0.066788i 0.350492 + 0.642638i 0.579629 + 0.048102i 0.600170 + 0.636938i ] * I;
1327 %!
1328 %! [Q,R] = qr (AA);
1329 %! [Q,R] = qrdelete (Q, R, 3);
1330 %! assert (norm (vec (Q'*Q - eye (5)), Inf) < 16*eps);
1331 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1332 %! assert (norm (vec (Q*R - [AA(:,1:2) AA(:,4)]), Inf) < norm (AA)*1e1*eps);
1333 %!
1334 %!test
1335 %! AA = [0.091364 0.613038 0.027504 0.999083;
1336 %! 0.594638 0.425302 0.562834 0.603537;
1337 %! 0.383594 0.291238 0.742073 0.085574;
1338 %! 0.265712 0.268003 0.783553 0.238409;
1339 %! 0.669966 0.743851 0.457255 0.445057 ];
1340 %!
1341 %! [Q,R] = qr (AA);
1342 %! [Q,R] = qrdelete (Q, R, 3, "row");
1343 %! assert (norm (vec (Q'*Q - eye (4)), Inf) < 1e1*eps);
1344 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1345 %! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf) < norm (AA)*1e1*eps);
1346 %!
1347 %!test
1348 %! AA = [0.364554 + 0.993117i 0.669818 + 0.510234i 0.426568 + 0.041337i 0.847051 + 0.233291i;
1349 %! 0.049600 + 0.242783i 0.448946 + 0.484022i 0.141155 + 0.074420i 0.446746 + 0.392706i;
1350 %! 0.581922 + 0.657416i 0.581460 + 0.030016i 0.219909 + 0.447288i 0.201144 + 0.069132i;
1351 %! 0.694986 + 0.000571i 0.682327 + 0.841712i 0.807537 + 0.166086i 0.192767 + 0.358098i;
1352 %! 0.945002 + 0.066788i 0.350492 + 0.642638i 0.579629 + 0.048102i 0.600170 + 0.636938i ] * I;
1353 %!
1354 %! [Q,R] = qr (AA);
1355 %! [Q,R] = qrdelete (Q, R, 3, "row");
1356 %! assert (norm (vec (Q'*Q - eye (4)), Inf) < 1e1*eps);
1357 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1358 %! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf) < norm (AA)*1e1*eps);
1359 
1360 %!test
1361 %! AA = single ([0.091364 0.613038 0.027504 0.999083;
1362 %! 0.594638 0.425302 0.562834 0.603537;
1363 %! 0.383594 0.291238 0.742073 0.085574;
1364 %! 0.265712 0.268003 0.783553 0.238409;
1365 %! 0.669966 0.743851 0.457255 0.445057 ]);
1366 %!
1367 %! [Q,R] = qr (AA);
1368 %! [Q,R] = qrdelete (Q, R, 3);
1369 %! assert (norm (vec (Q'*Q - eye (5,"single")), Inf) < 1e1*eps ("single"));
1370 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1371 %! assert (norm (vec (Q*R - [AA(:,1:2) AA(:,4)]), Inf) < norm (AA)*1e1*eps ("single"));
1372 %!
1373 %!test
1374 %! AA = single ([0.364554 + 0.993117i 0.669818 + 0.510234i 0.426568 + 0.041337i 0.847051 + 0.233291i;
1375 %! 0.049600 + 0.242783i 0.448946 + 0.484022i 0.141155 + 0.074420i 0.446746 + 0.392706i;
1376 %! 0.581922 + 0.657416i 0.581460 + 0.030016i 0.219909 + 0.447288i 0.201144 + 0.069132i;
1377 %! 0.694986 + 0.000571i 0.682327 + 0.841712i 0.807537 + 0.166086i 0.192767 + 0.358098i;
1378 %! 0.945002 + 0.066788i 0.350492 + 0.642638i 0.579629 + 0.048102i 0.600170 + 0.636938i ]) * I;
1379 %!
1380 %! [Q,R] = qr (AA);
1381 %! [Q,R] = qrdelete (Q, R, 3);
1382 %! assert (norm (vec (Q'*Q - eye (5,"single")), Inf) < 1e1*eps ("single"));
1383 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1384 %! assert (norm (vec (Q*R - [AA(:,1:2) AA(:,4)]), Inf) < norm (AA)*1e1*eps ("single"));
1385 %!
1386 %!test
1387 %! AA = single ([0.091364 0.613038 0.027504 0.999083;
1388 %! 0.594638 0.425302 0.562834 0.603537;
1389 %! 0.383594 0.291238 0.742073 0.085574;
1390 %! 0.265712 0.268003 0.783553 0.238409;
1391 %! 0.669966 0.743851 0.457255 0.445057 ]);
1392 %!
1393 %! [Q,R] = qr (AA);
1394 %! [Q,R] = qrdelete (Q, R, 3, "row");
1395 %! assert (norm (vec (Q'*Q - eye (4,"single")), Inf) < 1.5e1*eps ("single"));
1396 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1397 %! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf) < norm (AA)*1e1*eps ("single"));
1398 %!testif HAVE_QRUPDATE
1399 %! ## Same test as above but with more precicision
1400 %! AA = single ([0.091364 0.613038 0.027504 0.999083;
1401 %! 0.594638 0.425302 0.562834 0.603537;
1402 %! 0.383594 0.291238 0.742073 0.085574;
1403 %! 0.265712 0.268003 0.783553 0.238409;
1404 %! 0.669966 0.743851 0.457255 0.445057 ]);
1405 %!
1406 %! [Q,R] = qr (AA);
1407 %! [Q,R] = qrdelete (Q, R, 3, "row");
1408 %! assert (norm (vec (Q'*Q - eye (4,"single")), Inf) < 1e1*eps ("single"));
1409 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1410 %! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf) < norm (AA)*1e1*eps ("single"));
1411 %!
1412 %!test
1413 %! AA = single ([0.364554 + 0.993117i 0.669818 + 0.510234i 0.426568 + 0.041337i 0.847051 + 0.233291i;
1414 %! 0.049600 + 0.242783i 0.448946 + 0.484022i 0.141155 + 0.074420i 0.446746 + 0.392706i;
1415 %! 0.581922 + 0.657416i 0.581460 + 0.030016i 0.219909 + 0.447288i 0.201144 + 0.069132i;
1416 %! 0.694986 + 0.000571i 0.682327 + 0.841712i 0.807537 + 0.166086i 0.192767 + 0.358098i;
1417 %! 0.945002 + 0.066788i 0.350492 + 0.642638i 0.579629 + 0.048102i 0.600170 + 0.636938i ]) * I;
1418 %!
1419 %! [Q,R] = qr (AA);
1420 %! [Q,R] = qrdelete (Q, R, 3, "row");
1421 %! assert (norm (vec (Q'*Q - eye (4,"single")), Inf) < 1e1*eps ("single"));
1422 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1423 %! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf) < norm (AA)*1e1*eps ("single"));
1424 */
1425 
1426 DEFUN_DLD (qrshift, args, ,
1427  "-*- texinfo -*-\n\
1428 @deftypefn {Loadable Function} {[@var{Q1}, @var{R1}] =} qrshift (@var{Q}, @var{R}, @var{i}, @var{j})\n\
1429 Given a QR@tie{}factorization of a real or complex matrix\n\
1430 @w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and\n\
1431 @var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization\n\
1432 of @w{@var{A}(:,p)}, where @w{p} is the permutation @*\n\
1433 @code{p = [1:i-1, shift(i:j, 1), j+1:n]} if @w{@var{i} < @var{j}} @*\n\
1434  or @*\n\
1435 @code{p = [1:j-1, shift(j:i,-1), i+1:n]} if @w{@var{j} < @var{i}}. @*\n\
1436 \n\
1437 @seealso{qr, qrupdate, qrinsert, qrdelete}\n\
1438 @end deftypefn")
1439 {
1440  octave_idx_type nargin = args.length ();
1441  octave_value_list retval;
1442 
1443  if (nargin != 4)
1444  {
1445  print_usage ();
1446  return retval;
1447  }
1448 
1449  octave_value argq = args(0);
1450  octave_value argr = args(1);
1451  octave_value argi = args(2);
1452  octave_value argj = args(3);
1453 
1454  if (argq.is_numeric_type () && argr.is_numeric_type ())
1455  {
1456  if (check_qr_dims (argq, argr, true))
1457  {
1458  if (check_index (argi) && check_index (argj))
1459  {
1460  octave_idx_type i = argi.int_value ();
1461  octave_idx_type j = argj.int_value ();
1462 
1463  if (argq.is_real_type ()
1464  && argr.is_real_type ())
1465  {
1466  // all real case
1467  if (argq.is_single_type ()
1468  && argr.is_single_type ())
1469  {
1470  FloatMatrix Q = argq.float_matrix_value ();
1471  FloatMatrix R = argr.float_matrix_value ();
1472 
1473  FloatQR fact (Q, R);
1474  fact.shift_cols (i-1, j-1);
1475 
1476  retval(1) = get_qr_r (fact);
1477  retval(0) = fact.Q ();
1478  }
1479  else
1480  {
1481  Matrix Q = argq.matrix_value ();
1482  Matrix R = argr.matrix_value ();
1483 
1484  QR fact (Q, R);
1485  fact.shift_cols (i-1, j-1);
1486 
1487  retval(1) = get_qr_r (fact);
1488  retval(0) = fact.Q ();
1489  }
1490  }
1491  else
1492  {
1493  // complex case
1494  if (argq.is_single_type ()
1495  && argr.is_single_type ())
1496  {
1499 
1500  FloatComplexQR fact (Q, R);
1501  fact.shift_cols (i-1, j-1);
1502 
1503  retval(1) = get_qr_r (fact);
1504  retval(0) = fact.Q ();
1505  }
1506  else
1507  {
1509  ComplexMatrix R = argr.complex_matrix_value ();
1510 
1511  ComplexQR fact (Q, R);
1512  fact.shift_cols (i-1, j-1);
1513 
1514  retval(1) = get_qr_r (fact);
1515  retval(0) = fact.Q ();
1516  }
1517  }
1518  }
1519  else
1520  error ("qrshift: invalid index I or J");
1521  }
1522  else
1523  error ("qrshift: dimensions mismatch");
1524  }
1525  else
1526  error ("qrshift: Q and R must be numeric");
1527 
1528  return retval;
1529 }
1530 
1531 /*
1532 %!test
1533 %! AA = A.';
1534 %! i = 2; j = 4; p = [1:i-1, shift(i:j,-1), j+1:5];
1535 %!
1536 %! [Q,R] = qr (AA);
1537 %! [Q,R] = qrshift (Q, R, i, j);
1538 %! assert (norm (vec (Q'*Q - eye (3)), Inf) < 1e1*eps);
1539 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1540 %! assert (norm (vec (Q*R - AA(:,p)), Inf) < norm (AA)*1e1*eps);
1541 %!
1542 %! j = 2; i = 4; p = [1:j-1, shift(j:i,+1), i+1:5];
1543 %!
1544 %! [Q,R] = qr (AA);
1545 %! [Q,R] = qrshift (Q, R, i, j);
1546 %! assert (norm (vec (Q'*Q - eye (3)), Inf) < 1e1*eps);
1547 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1548 %! assert (norm (vec (Q*R - AA(:,p)), Inf) < norm (AA)*1e1*eps);
1549 %!
1550 %!test
1551 %! AA = Ac.';
1552 %! i = 2; j = 4; p = [1:i-1, shift(i:j,-1), j+1:5];
1553 %!
1554 %! [Q,R] = qr (AA);
1555 %! [Q,R] = qrshift (Q, R, i, j);
1556 %! assert (norm (vec (Q'*Q - eye (3)), Inf) < 1e1*eps);
1557 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1558 %! assert (norm (vec (Q*R - AA(:,p)), Inf) < norm (AA)*1e1*eps);
1559 %!
1560 %! j = 2; i = 4; p = [1:j-1, shift(j:i,+1), i+1:5];
1561 %!
1562 %! [Q,R] = qr (AA);
1563 %! [Q,R] = qrshift (Q, R, i, j);
1564 %! assert (norm (vec (Q'*Q - eye (3)), Inf) < 1e1*eps);
1565 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1566 %! assert (norm (vec (Q*R - AA(:,p)), Inf) < norm (AA)*1e1*eps);
1567 
1568 %!test
1569 %! AA = single (A).';
1570 %! i = 2; j = 4; p = [1:i-1, shift(i:j,-1), j+1:5];
1571 %!
1572 %! [Q,R] = qr (AA);
1573 %! [Q,R] = qrshift (Q, R, i, j);
1574 %! assert (norm (vec (Q'*Q - eye (3,"single")), Inf) < 1e1*eps ("single"));
1575 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1576 %! assert (norm (vec (Q*R - AA(:,p)), Inf) < norm (AA)*1e1*eps ("single"));
1577 %!
1578 %! j = 2; i = 4; p = [1:j-1, shift(j:i,+1), i+1:5];
1579 %!
1580 %! [Q,R] = qr (AA);
1581 %! [Q,R] = qrshift (Q, R, i, j);
1582 %! assert (norm (vec (Q'*Q - eye (3,"single")), Inf) < 1e1*eps ("single"));
1583 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1584 %! assert (norm (vec (Q*R - AA(:,p)), Inf) < norm (AA)*1e1*eps ("single"));
1585 %!
1586 %!test
1587 %! AA = single (Ac).';
1588 %! i = 2; j = 4; p = [1:i-1, shift(i:j,-1), j+1:5];
1589 %!
1590 %! [Q,R] = qr (AA);
1591 %! [Q,R] = qrshift (Q, R, i, j);
1592 %! assert (norm (vec (Q'*Q - eye (3,"single")), Inf) < 1e1*eps ("single"));
1593 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1594 %! assert (norm (vec (Q*R - AA(:,p)), Inf) < norm (AA)*1e1*eps ("single"));
1595 %!
1596 %! j = 2; i = 4; p = [1:j-1, shift(j:i,+1), i+1:5];
1597 %!
1598 %! [Q,R] = qr (AA);
1599 %! [Q,R] = qrshift (Q, R, i, j);
1600 %! assert (norm (vec (Q'*Q - eye (3,"single")), Inf) < 1e1*eps ("single"));
1601 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1602 %! assert (norm (vec (Q*R - AA(:,p)), Inf) < norm (AA)*1e1*eps ("single"));
1603 */
void insert_col(const ColumnVector &u, octave_idx_type j)
Definition: dbleQR.cc:254
void delete_row(octave_idx_type j)
Definition: dbleQR.cc:424
static const type std
Definition: dbleQR.h:43
bool is_real_type(void) const
Definition: ov.h:651
void gripe_wrong_type_arg(const char *name, const char *s, bool is_error)
Definition: gripes.cc:135
int ndims(void) const
Definition: ov.h:479
octave_idx_type rows(void) const
Definition: ov.h:473
qr_type Q(void) const
Definition: base-qr.h:63
FloatComplexMatrix float_complex_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:795
OCTINTERP_API void print_usage(void)
Definition: defun.cc:51
static octave_value get_qr_r(const base_qr< MT > &fact)
Definition: qr.cc:49
octave_idx_type length(void) const
Definition: oct-obj.h:89
void shift_cols(octave_idx_type i, octave_idx_type j)
Definition: CmplxQR.cc:449
bool is_scalar_type(void) const
Definition: ov.h:657
bool is_numeric_type(void) const
Definition: ov.h:663
void update(const ComplexColumnVector &u, const ComplexColumnVector &v)
Definition: CmplxQR.cc:209
int int_value(bool req_int=false, bool frc_str_conv=false) const
Definition: ov.h:730
RowVector Pvec(void) const
Definition: CmplxQRP.cc:104
void error(const char *fmt,...)
Definition: error.cc:476
void delete_row(octave_idx_type j)
Definition: CmplxQR.cc:426
RowVector row(octave_idx_type i) const
Definition: dMatrix.cc:639
void insert_row(const FloatRowVector &u, octave_idx_type j)
Definition: floatQR.cc:398
PermMatrix P(void) const
Definition: floatQRP.h:59
int empty_arg(const char *, octave_idx_type nr, octave_idx_type nc)
Definition: utils.cc:246
ComplexRowVector row(octave_idx_type i) const
Definition: CMatrix.cc:990
static const type raw
Definition: dbleQR.h:44
void shift_cols(octave_idx_type i, octave_idx_type j)
Definition: dbleQR.cc:446
FloatRowVector row(octave_idx_type i) const
Definition: fMatrix.cc:646
void shift_cols(octave_idx_type i, octave_idx_type j)
Definition: floatQR.cc:444
bool regular(void) const
Definition: base-qr.cc:65
void insert_row(const RowVector &u, octave_idx_type j)
Definition: dbleQR.cc:400
void delete_row(octave_idx_type j)
Definition: floatQR.cc:422
octave_idx_type columns(void) const
Definition: ov.h:475
FloatRowVector Pvec(void) const
Definition: floatQRP.cc:101
static const type economy
Definition: dbleQR.h:45
bool is_sparse_type(void) const
Definition: ov.h:666
RowVector Pvec(void) const
Definition: dbleQRP.cc:101
int error_state
Definition: error.cc:101
Array< octave_idx_type > octave_idx_type_vector_value(bool req_int=false, bool frc_str_conv=false, bool frc_vec_conv=false) const
Definition: ov.cc:1792
bool is_complex_type(void) const
Definition: ov.h:654
void shift_cols(octave_idx_type i, octave_idx_type j)
Definition: fCmplxQR.cc:455
qr_type_t
Definition: base-qr.h:30
void insert_row(const ComplexRowVector &u, octave_idx_type j)
Definition: CmplxQR.cc:402
void delete_col(octave_idx_type j)
Definition: fCmplxQR.cc:339
void insert_col(const FloatComplexColumnVector &u, octave_idx_type j)
Definition: fCmplxQR.cc:261
Definition: dMatrix.h:35
SparseComplexMatrix sparse_complex_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:824
Matrix matrix_value(bool frc_str_conv=false) const
Definition: ov.h:773
static bool check_index(const octave_value &i, bool vector_allowed=false)
Definition: qr.cc:757
Definition: dbleQRP.h:32
double arg(double x)
Definition: lo-mappers.h:37
FloatRowVector Pvec(void) const
Definition: fCmplxQRP.cc:104
void delete_col(octave_idx_type j)
Definition: floatQR.cc:329
void warning(const char *fmt,...)
Definition: error.cc:681
void insert_col(const FloatColumnVector &u, octave_idx_type j)
Definition: floatQR.cc:252
void update(const FloatColumnVector &u, const FloatColumnVector &v)
Definition: floatQR.cc:208
Definition: dbleQR.h:34
void insert_col(const ComplexColumnVector &u, octave_idx_type j)
Definition: CmplxQR.cc:256
qr_type R(void) const
Definition: base-qr.h:65
ComplexMatrix complex_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:791
void delete_col(octave_idx_type j)
Definition: CmplxQR.cc:333
FloatMatrix float_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:776
void delete_row(octave_idx_type j)
Definition: fCmplxQR.cc:432
void insert_row(const FloatComplexRowVector &u, octave_idx_type j)
Definition: fCmplxQR.cc:408
F77_RET_T const octave_idx_type const octave_idx_type const octave_idx_type double const octave_idx_type double const octave_idx_type double * Q
Definition: qz.cc:114
#define DEFUN_DLD(name, args_name, nargout_name, doc)
Definition: defun-dld.h:59
SparseMatrix sparse_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:820
bool is_single_type(void) const
Definition: ov.h:611
FloatComplexRowVector row(octave_idx_type i) const
Definition: fCMatrix.cc:996
void delete_col(octave_idx_type j)
Definition: dbleQR.cc:331
PermMatrix P(void) const
Definition: CmplxQRP.h:58
void update(const ColumnVector &u, const ColumnVector &v)
Definition: dbleQR.cc:210
F77_RET_T const double * x
PermMatrix P(void) const
Definition: dbleQRP.h:59
void update(const FloatComplexColumnVector &u, const FloatComplexColumnVector &v)
Definition: fCmplxQR.cc:212
static bool check_qr_dims(const octave_value &q, const octave_value &r, bool allow_ecf=false)
Definition: qr.cc:746
PermMatrix P(void) const
Definition: fCmplxQRP.h:58
bool is_integer_type(void) const
Definition: ov.h:648