GNU Octave  8.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
qr.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 //
3 // Copyright (C) 1996-2023 The Octave Project Developers
4 //
5 // See the file COPYRIGHT.md in the top-level directory of this
6 // distribution or <https://octave.org/copyright/>.
7 //
8 // This file is part of Octave.
9 //
10 // Octave is free software: you can redistribute it and/or modify it
11 // under the terms of the GNU General Public License as published by
12 // the Free Software Foundation, either version 3 of the License, or
13 // (at your option) any later version.
14 //
15 // Octave is distributed in the hope that it will be useful, but
16 // WITHOUT ANY WARRANTY; without even the implied warranty of
17 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 // GNU General Public License for more details.
19 //
20 // You should have received a copy of the GNU General Public License
21 // along with Octave; see the file COPYING. If not, see
22 // <https://www.gnu.org/licenses/>.
23 //
24 ////////////////////////////////////////////////////////////////////////
25 
26 #if defined (HAVE_CONFIG_H)
27 # include "config.h"
28 #endif
29 
30 #include <string>
31 
32 #include "MArray.h"
33 #include "Matrix.h"
34 #include "qr.h"
35 #include "qrp.h"
36 #include "sparse-qr.h"
37 
38 #include "defun.h"
39 #include "error.h"
40 #include "errwarn.h"
41 #include "ov.h"
42 #include "ovl.h"
43 
45 
46 /*
47 ## Restore all rand* "state" values
48 %!function restore_rand_states (state)
49 %! rand ("state", state.rand);
50 %! randn ("state", state.randn);
51 %!endfunction
52 
53 %!shared old_state, restore_state
54 %! ## Save and restore the states of both random number generators that are
55 %! ## tested by the unit tests in this file.
56 %! old_state.rand = rand ("state");
57 %! old_state.randn = randn ("state");
58 %! restore_state = onCleanup (@() restore_rand_states (old_state));
59 */
60 
61 template <typename MT>
62 static octave_value
63 get_qr_r (const math::qr<MT>& fact)
64 {
65  MT R = fact.R ();
66  if (R.issquare () && fact.regular ())
68  else
69  return R;
70 }
71 
72 template <typename T>
73 static typename math::qr<T>::type
74 qr_type (int nargout, bool economy)
75 {
76  if (nargout == 0 || nargout == 1)
77  return math::qr<T>::raw;
78  else if (economy)
79  return math::qr<T>::economy;
80  else
81  return math::qr<T>::std;
82 }
83 
84 // dense X
85 //
86 // [Q, R] = qr (X): form Q unitary and R upper triangular such
87 // that Q * R = X
88 //
89 // [Q, R] = qr (X, 0): form the economy decomposition such that if X is
90 // m by n then only the first n columns of Q are
91 // computed.
92 //
93 // [Q, R, P] = qr (X): form QRP factorization of X where
94 // P is a permutation matrix such that
95 // A * P = Q * R
96 //
97 // [Q, R, P] = qr (X, 0): form the economy decomposition with
98 // permutation vector P such that Q * R = X(:, P)
99 //
100 // qr (X) alone returns the output of the LAPACK routine dgeqrf, such
101 // that R = triu (qr (X))
102 //
103 // sparse X
104 //
105 // X = qr (A, B): if M < N, X is the minimum 2-norm solution of
106 // A\B. If M >= N, X is the least squares
107 // approximation of A\B. X is calculated by
108 // SPQR-function SuiteSparseQR_min2norm.
109 //
110 DEFUN (qr, args, nargout,
111  doc: /* -*- texinfo -*-
112 @deftypefn {} {[@var{Q}, @var{R}] =} qr (@var{A})
113 @deftypefnx {} {[@var{Q}, @var{R}, @var{P}] =} qr (@var{A})
114 @deftypefnx {} {@var{X} =} qr (@var{A}) # non-sparse A
115 @deftypefnx {} {@var{R} =} qr (@var{A}) # sparse A
116 @deftypefnx {} {@var{X} =} qr (@var{A}, @var{B}) # sparse A
117 @deftypefnx {} {[@var{C}, @var{R}] =} qr (@var{A}, @var{B})
118 @deftypefnx {} {[@dots{}] =} qr (@dots{}, 0)
119 @deftypefnx {} {[@dots{}] =} qr (@dots{}, "econ")
120 @deftypefnx {} {[@dots{}] =} qr (@dots{}, "vector")
121 @deftypefnx {} {[@dots{}] =} qr (@dots{}, "matrix")
122 @cindex QR factorization
123 Compute the QR@tie{}factorization of @var{A}, using standard @sc{lapack}
124 subroutines.
125 
126 The QR@tie{}factorization is
127 @tex
128 $QR = A$ where $Q$ is an orthogonal matrix and $R$ is upper triangular.
129 @end tex
130 @ifnottex
131 
132 @example
133 @var{Q} * @var{R} = @var{A}
134 @end example
135 
136 @noindent
137 where @var{Q} is an orthogonal matrix and @var{R} is upper triangular.
138 @end ifnottex
139 
140 For example, given the matrix @code{@var{A} = [1, 2; 3, 4]},
141 
142 @example
143 [@var{Q}, @var{R}] = qr (@var{A})
144 @end example
145 
146 @noindent
147 returns
148 
149 @example
150 @group
151 @var{Q} =
152 
153  -0.31623 -0.94868
154  -0.94868 0.31623
155 
156 @var{R} =
157 
158  -3.16228 -4.42719
159  0.00000 -0.63246
160 @end group
161 @end example
162 
163 @noindent
164 which multiplied together return the original matrix
165 
166 @example
167 @group
168 @var{Q} * @var{R}
169  @result{}
170  1.0000 2.0000
171  3.0000 4.0000
172 @end group
173 @end example
174 
175 If just a single return value is requested then it is either @var{R}, if
176 @var{A} is sparse, or @var{X}, such that @code{@var{R} = triu (@var{X})} if
177 @var{A} is full. (Note: unlike most commands, the single return value is not
178 the first return value when multiple values are requested.)
179 
180 If a third output @var{P} is requested, then @code{qr} calculates the permuted
181 QR@tie{}factorization
182 @tex
183 $QR = AP$ where $Q$ is an orthogonal matrix, $R$ is upper triangular, and $P$
184 is a permutation matrix.
185 @end tex
186 @ifnottex
187 
188 @example
189 @var{Q} * @var{R} = @var{A} * @var{P}
190 @end example
191 
192 @noindent
193 where @var{Q} is an orthogonal matrix, @var{R} is upper triangular, and
194 @var{P} is a permutation matrix.
195 @end ifnottex
196 
197 If @var{A} is dense, the permuted QR@tie{}factorization has the additional
198 property that the diagonal entries of @var{R} are ordered by decreasing
199 magnitude. In other words, @code{abs (diag (@var{R}))} will be ordered
200 from largest to smallest.
201 
202 If @var{A} is sparse, @var{P} is a fill-reducing ordering of the columns
203 of @var{A}. In that case, the diagonal entries of @var{R} are not ordered by
204 decreasing magnitude.
205 
206 For example, given the matrix @code{@var{A} = [1, 2; 3, 4]},
207 
208 @example
209 [@var{Q}, @var{R}, @var{P}] = qr (@var{A})
210 @end example
211 
212 @noindent
213 returns
214 
215 @example
216 @group
217 @var{Q} =
218 
219  -0.44721 -0.89443
220  -0.89443 0.44721
221 
222 @var{R} =
223 
224  -4.47214 -3.13050
225  0.00000 0.44721
226 
227 @var{P} =
228 
229  0 1
230  1 0
231 @end group
232 @end example
233 
234 If the input matrix @var{A} is sparse, the sparse QR@tie{}factorization
235 is computed by using @sc{SPQR} or @sc{cxsparse} (e.g., if @sc{SPQR} is not
236 available). Because the matrix @var{Q} is, in general, a full matrix, it is
237 recommended to request only one return value @var{R}. In that case, the
238 computation avoids the construction of @var{Q} and returns a sparse @var{R}
239 such that @code{@var{R} = chol (@var{A}' * @var{A})}.
240 
241 If @var{A} is dense, an additional matrix @var{B} is supplied and two
242 return values are requested, then @code{qr} returns @var{C}, where
243 @code{@var{C} = @var{Q}' * @var{B}}. This allows the least squares
244 approximation of @code{@var{A} \ @var{B}} to be calculated as
245 
246 @example
247 @group
248 [@var{C}, @var{R}] = qr (@var{A}, @var{B})
249 @var{X} = @var{R} \ @var{C}
250 @end group
251 @end example
252 
253 If @var{A} is a sparse @nospell{MxN} matrix and an additional matrix @var{B} is
254 supplied, one or two return values are possible. If one return value @var{X}
255 is requested and M < N, then @var{X} is the minimum 2-norm solution of
256 @w{@code{@var{A} \ @var{B}}}. If M >= N, @var{X} is the least squares
257 approximation @w{of @code{@var{A} \ @var{B}}}. If two return values are
258 requested, @var{C} and @var{R} have the same meaning as in the dense case
259 (@var{C} is dense and @var{R} is sparse). The version with one return
260 parameter should be preferred because it uses less memory and can handle
261 rank-deficient matrices better.
262 
263 If the final argument is the string @qcode{"vector"} then @var{P} is a
264 permutation vector (of the columns of @var{A}) instead of a permutation
265 matrix. In this case, the defining relationship is:
266 
267 @example
268 @var{Q} * @var{R} = @var{A}(:, @var{P})
269 @end example
270 
271 The default, however, is to return a permutation matrix and this may be
272 explicitly specified by using a final argument of @qcode{"matrix"}.
273 
274 If the final argument is the scalar 0 or the string @qcode{"econ"}, an economy
275 factorization is returned. If the original matrix @var{A} has size
276 @nospell{MxN} and M > N, then the economy factorization will calculate just N
277 rows in @var{R} and N columns in @var{Q} and omit the zeros in @var{R}. If M
278 @leq{} N, there is no difference between the economy and standard
279 factorizations. When calculating an economy factorization and @var{A} is
280 dense, the output @var{P} is always a vector rather than a matrix. If @var{A}
281 is sparse, output @var{P} is a sparse permutation matrix.
282 
283 Background: The QR factorization has applications in the solution of least
284 squares problems
285 @tex
286 $$
287 \min_x \left\Vert A x - b \right\Vert_2
288 $$
289 @end tex
290 @ifnottex
291 
292 @example
293 min norm (A*x - b)
294 @end example
295 
296 @end ifnottex
297 for overdetermined systems of equations (i.e.,
298 @tex
299 $A$
300 @end tex
301 @ifnottex
302 @var{A}
303 @end ifnottex
304 is a tall, thin matrix).
305 
306 The permuted QR@tie{}factorization
307 @code{[@var{Q}, @var{R}, @var{P}] = qr (@var{A})} allows the construction of an
308 orthogonal basis of @code{span (A)}.
309 
310 @seealso{chol, hess, lu, qz, schur, svd, qrupdate, qrinsert, qrdelete, qrshift}
311 @end deftypefn */)
312 {
313  int nargin = args.length ();
314 
315  if (nargin < 1 || nargin > 3)
316  print_usage ();
317 
318  octave_value_list retval;
319 
320  octave_value arg = args(0);
321 
322  bool economy = false;
323  bool is_cmplx = false;
324  bool have_b = false;
325  bool vector_p = false;
326 
327  if (arg.iscomplex ())
328  is_cmplx = true;
329  if (nargin > 1)
330  {
331  have_b = true;
332  if (args(nargin-1).is_scalar_type ())
333  {
334  int val = args(nargin-1).int_value ();
335  if (val == 0)
336  {
337  economy = true;
338  have_b = (nargin > 2);
339  }
340  else if (nargin == 3) // argument 3 should be 0 or a string
341  print_usage ();
342  }
343  else if (args(nargin-1).is_string ())
344  {
345  std::string str = args(nargin-1).string_value ();
346  if (str == "vector")
347  vector_p = true;
348  else if (str == "econ")
349  {
350  economy = true;
351  have_b = (nargin > 2);
352  }
353  else if (str != "matrix")
354  error ("qr: option string must be 'econ' or 'matrix' or " \
355  "'vector', not \"%s\"", str.c_str ());
356  have_b = (nargin > 2);
357  }
358  else if (! args(nargin-1).is_matrix_type ())
359  err_wrong_type_arg ("qr", args(nargin-1));
360  else if (nargin == 3) // should be caught by is_scalar_type or is_string
361  print_usage ();
362 
363  if (have_b && args(1).iscomplex ())
364  is_cmplx = true;
365  }
366 
367  if (arg.issparse ())
368  {
369  if (nargout > 3)
370  error ("qr: too many output arguments");
371 
372  if (is_cmplx)
373  {
374  if (have_b && nargout == 1)
375  {
376  octave_idx_type info;
377 
378  if (! args(1).issparse () && args(1).iscomplex ())
379  retval = ovl
380  (math::sparse_qr<SparseComplexMatrix>::solve
383  args(1).complex_matrix_value (), info));
384  else if (args(1).issparse () && args(1).iscomplex ())
385  retval = ovl
386  (math::sparse_qr<SparseComplexMatrix>::solve
389  args(1).sparse_complex_matrix_value (), info));
390  else if (! args(1).issparse () && ! args(1).iscomplex ())
391  retval = ovl
392  (math::sparse_qr<SparseComplexMatrix>::solve
395  args(1).matrix_value (), info));
396  else if (args(1).issparse () && ! args(1).iscomplex ())
397  retval = ovl
398  (math::sparse_qr<SparseComplexMatrix>::solve
401  args(1).sparse_matrix_value (), info));
402  else
403  error ("qr: b is not valid");
404  }
405  else if (have_b && nargout == 2)
406  {
408  q (arg.sparse_complex_matrix_value (), 0);
409  retval = ovl (q.C (args(1).complex_matrix_value (), economy),
410  q.R (economy));
411  }
412  else if (have_b && nargout == 3)
413  {
415  q (arg.sparse_complex_matrix_value ());
416  if (vector_p)
417  retval = ovl (q.C (args(1).complex_matrix_value (), economy),
418  q.R (economy), q.E ());
419  else
420  retval = ovl (q.C (args(1).complex_matrix_value (), economy),
421  q.R (economy), q.E_MAT ());
422  }
423  else
424  {
425  if (nargout > 2)
426  {
428  q (arg.sparse_complex_matrix_value ());
429  if (vector_p)
430  retval = ovl (q.Q (economy), q.R (economy), q.E ());
431  else
432  retval = ovl (q.Q (economy), q.R (economy),
433  q.E_MAT ());
434  }
435  else if (nargout > 1)
436  {
438  q (arg.sparse_complex_matrix_value (), 0);
439  retval = ovl (q.Q (economy), q.R (economy));
440  }
441  else
442  {
444  q (arg.sparse_complex_matrix_value (), 0);
445  retval = ovl (q.R (economy));
446  }
447  }
448  }
449  else
450  {
451  if (have_b && nargout == 1)
452  {
453  octave_idx_type info;
454  if (args(1).issparse () && ! args(1).iscomplex ())
455  retval = ovl (math::sparse_qr<SparseMatrix>::solve
457  (arg.sparse_matrix_value (),
458  args (1).sparse_matrix_value (), info));
459  else if (! args(1).issparse () && args(1).iscomplex ())
460  retval = ovl (math::sparse_qr<SparseMatrix>::solve
462  (arg.sparse_matrix_value (),
463  args (1).complex_matrix_value (), info));
464  else if (! args(1).issparse () && ! args(1).iscomplex ())
465  retval = ovl (math::sparse_qr<SparseMatrix>::solve
467  (arg.sparse_matrix_value (),
468  args (1).matrix_value (), info));
469  else if (args(1).issparse () && args(1).iscomplex ())
470  retval = ovl (math::sparse_qr<SparseMatrix>::solve
472  (arg.sparse_matrix_value (),
473  args(1).sparse_complex_matrix_value (),
474  info));
475  else
476  error ("qr: b is not valid");
477  }
478  else if (have_b && nargout == 2)
479  {
481  q (arg.sparse_matrix_value (), 0);
482  retval = ovl (q.C (args(1).matrix_value (), economy),
483  q.R (economy));
484  }
485  else if (have_b && nargout == 3)
486  {
488  q (arg.sparse_matrix_value ());
489  if (vector_p)
490  retval = ovl (q.C (args(1).matrix_value (), economy),
491  q.R (economy), q.E ());
492  else
493  retval = ovl (q.C (args(1).matrix_value (), economy),
494  q.R (economy), q.E_MAT ());
495  }
496 
497  else
498  {
499  if (nargout > 2)
500  {
502  q (arg.sparse_matrix_value ());
503  if (vector_p)
504  retval = ovl (q.Q (economy), q.R (economy), q.E ());
505  else
506  retval = ovl (q.Q (economy), q.R (economy),
507  q.E_MAT ());
508  }
509  else if (nargout > 1)
510  {
512  q (arg.sparse_matrix_value (), 0);
513  retval = ovl (q.Q (economy), q.R (economy));
514  }
515  else
516  {
518  q (arg.sparse_matrix_value (), 0);
519  retval = ovl (q.R (economy));
520  }
521  }
522  }
523  }
524  else
525  {
526  if (have_b && nargout > 2)
527  error ("qr: too many output arguments for dense A with B");
528 
529  if (arg.is_single_type ())
530  {
531  if (arg.isreal ())
532  {
533  math::qr<FloatMatrix>::type type
534  = qr_type<FloatMatrix> (nargout, economy);
535 
537 
538  switch (nargout)
539  {
540  case 0:
541  case 1:
542  {
543  math::qr<FloatMatrix> fact (m, type);
544  retval = ovl (fact.R ());
545  }
546  break;
547 
548  case 2:
549  {
550  math::qr<FloatMatrix> fact (m, type);
551  retval = ovl (fact.Q (), get_qr_r (fact));
552  if (have_b)
553  {
554  if (is_cmplx)
555  retval(0) = fact.Q ().transpose ()
556  * args(1).float_complex_matrix_value ();
557  else
558  retval(0) = fact.Q ().transpose ()
559  * args(1).float_matrix_value ();
560  }
561  }
562  break;
563 
564  default:
565  {
566  math::qrp<FloatMatrix> fact (m, type);
567 
568  if (economy || vector_p)
569  retval = ovl (fact.Q (), get_qr_r (fact), fact.Pvec ());
570  else
571  retval = ovl (fact.Q (), get_qr_r (fact), fact.P ());
572  }
573  break;
574  }
575  }
576  else if (arg.iscomplex ())
577  {
578  math::qr<FloatComplexMatrix>::type type
579  = qr_type<FloatComplexMatrix> (nargout, economy);
580 
582 
583  switch (nargout)
584  {
585  case 0:
586  case 1:
587  {
588  math::qr<FloatComplexMatrix> fact (m, type);
589  retval = ovl (fact.R ());
590  }
591  break;
592 
593  case 2:
594  {
595  math::qr<FloatComplexMatrix> fact (m, type);
596  retval = ovl (fact.Q (), get_qr_r (fact));
597  if (have_b)
598  retval(0) = conj (fact.Q ().transpose ())
599  * args(1).float_complex_matrix_value ();
600  }
601  break;
602 
603  default:
604  {
605  math::qrp<FloatComplexMatrix> fact (m, type);
606  if (economy || vector_p)
607  retval = ovl (fact.Q (), get_qr_r (fact), fact.Pvec ());
608  else
609  retval = ovl (fact.Q (), get_qr_r (fact), fact.P ());
610  }
611  break;
612  }
613  }
614  }
615  else
616  {
617  if (arg.isreal ())
618  {
619  math::qr<Matrix>::type type
620  = qr_type<Matrix> (nargout, economy);
621 
622  Matrix m = arg.matrix_value ();
623 
624  switch (nargout)
625  {
626  case 0:
627  case 1:
628  {
629  math::qr<Matrix> fact (m, type);
630  retval = ovl (fact.R ());
631  }
632  break;
633 
634  case 2:
635  {
636  math::qr<Matrix> fact (m, type);
637  retval = ovl (fact.Q (), get_qr_r (fact));
638  if (have_b)
639  {
640  if (is_cmplx)
641  retval(0) = fact.Q ().transpose ()
642  * args(1).complex_matrix_value ();
643  else
644  retval(0) = fact.Q ().transpose ()
645  * args(1).matrix_value ();
646  }
647  }
648  break;
649 
650  default:
651  {
652  math::qrp<Matrix> fact (m, type);
653  if (economy || vector_p)
654  retval = ovl (fact.Q (), get_qr_r (fact), fact.Pvec ());
655  else
656  retval = ovl (fact.Q (), get_qr_r (fact), fact.P ());
657  }
658  break;
659  }
660  }
661  else if (arg.iscomplex ())
662  {
663  math::qr<ComplexMatrix>::type type
664  = qr_type<ComplexMatrix> (nargout, economy);
665 
667 
668  switch (nargout)
669  {
670  case 0:
671  case 1:
672  {
673  math::qr<ComplexMatrix> fact (m, type);
674  retval = ovl (fact.R ());
675  }
676  break;
677 
678  case 2:
679  {
680  math::qr<ComplexMatrix> fact (m, type);
681  retval = ovl (fact.Q (), get_qr_r (fact));
682  if (have_b)
683  retval(0) = conj (fact.Q ().transpose ())
684  * args(1).complex_matrix_value ();
685  }
686  break;
687 
688  default:
689  {
690  math::qrp<ComplexMatrix> fact (m, type);
691  if (economy || vector_p)
692  retval = ovl (fact.Q (), get_qr_r (fact), fact.Pvec ());
693  else
694  retval = ovl (fact.Q (), get_qr_r (fact), fact.P ());
695  }
696  break;
697  }
698  }
699  else
700  err_wrong_type_arg ("qr", arg);
701  }
702  }
703 
704  return retval;
705 }
706 
707 /*
708 %!test
709 %! a = [0, 2, 1; 2, 1, 2];
710 %!
711 %! [q, r] = qr (a);
712 %! [qe, re] = qr (a, 0);
713 %! [qe2, re2] = qr (a, "econ");
714 %!
715 %! assert (q * r, a, sqrt (eps));
716 %! assert (qe * re, a, sqrt (eps));
717 %! assert (qe2 * re2, a, sqrt (eps));
718 
719 %!test
720 %! a = [0, 2, 1; 2, 1, 2];
721 %!
722 %! [q, r] = qr (a);
723 %! [qe, re] = qr (a, 0);
724 %! [qe2, re2] = qr (a, "econ");
725 %!
726 %! assert (q * r, a, sqrt (eps));
727 %! assert (qe * re, a, sqrt (eps));
728 %! assert (qe2 * re2, a, sqrt (eps));
729 
730 %!test
731 %! a = [0, 2, 1; 2, 1, 2];
732 %!
733 %! [q, r, p] = qr (a); # FIXME: not giving right dimensions.
734 %! [qe, re, pe] = qr (a, 0);
735 %! [qe2, re2, pe2] = qr (a, "econ");
736 %!
737 %! assert (q * r, a * p, sqrt (eps));
738 %! assert (qe * re, a(:, pe), sqrt (eps));
739 %! assert (qe2 * re2, a(:, pe2), sqrt (eps));
740 
741 %!test
742 %! a = [0, 2; 2, 1; 1, 2];
743 %!
744 %! [q, r] = qr (a);
745 %! [qe, re] = qr (a, 0);
746 %! [qe2, re2] = qr (a, "econ");
747 %!
748 %! assert (q * r, a, sqrt (eps));
749 %! assert (qe * re, a, sqrt (eps));
750 %! assert (qe2 * re2, a, sqrt (eps));
751 
752 %!test
753 %! a = [0, 2; 2, 1; 1, 2];
754 %!
755 %! [q, r, p] = qr (a);
756 %! [qe, re, pe] = qr (a, 0);
757 %! [qe2, re2, pe2] = qr (a, "econ");
758 %!
759 %! assert (q * r, a * p, sqrt (eps));
760 %! assert (qe * re, a(:, pe), sqrt (eps));
761 %! assert (qe2 * re2, a(:, pe2), sqrt (eps));
762 
763 %!test
764 %! a = [0, 2, 1; 2, 1, 2; 3, 1, 2];
765 %! b = [1, 3, 2; 1, 1, 0; 3, 0, 2];
766 %!
767 %! [q, r] = qr (a);
768 %! [c, re] = qr (a, b);
769 %!
770 %! assert (r, re, sqrt (eps));
771 %! assert (q'*b, c, sqrt (eps));
772 
773 %!test
774 %! a = [0, 2, i; 2, 1, 2; 3, 1, 2];
775 %! b = [1, 3, 2; 1, i, 0; 3, 0, 2];
776 %!
777 %! [q, r] = qr (a);
778 %! [c, re] = qr (a, b);
779 %!
780 %! assert (r, re, sqrt (eps));
781 %! assert (q'*b, c, sqrt (eps));
782 
783 %!test
784 %! a = [0, 2, i; 2, 1, 2; 3, 1, 2];
785 %! b = [1, 3, 2; 1, 1, 0; 3, 0, 2];
786 %!
787 %! [q, r] = qr (a);
788 %! [c, re] = qr (a, b);
789 %!
790 %! assert (r, re, sqrt (eps));
791 %! assert (q'*b, c, sqrt (eps));
792 
793 %!test
794 %! a = [0, 2, 1; 2, 1, 2; 3, 1, 2];
795 %! b = [1, 3, 2; 1, i, 0; 3, 0, 2];
796 %!
797 %! [q, r] = qr (a);
798 %! [c, re] = qr (a, b);
799 %!
800 %! assert (r, re, sqrt (eps));
801 %! assert (q'*b, c, sqrt (eps));
802 
803 ## Empty matrices
804 %!test
805 %! assert (qr (zeros (0, 0)), zeros (0, 0))
806 %! assert (qr (zeros (1, 0)), zeros (1, 0))
807 %! assert (qr (zeros (0, 1)), zeros (0, 1))
808 
809 %!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE") <*63069>
810 %! assert (qr (sparse (0, 0)), sparse (0, 0))
811 %! assert (qr (sparse (1, 0)), sparse (1, 0))
812 %! assert (qr (sparse (0, 1)), sparse (0, 1))
813 
814 %!error qr ()
815 %!error qr ([1, 2; 3, 4], 0, 2)
816 %!error <option string must be .*, not "foo"> qr (magic (3), "foo")
817 %!error <option string must be .*, not "foo"> qr (magic (3), rand (3, 1), "foo")
818 %!error <too many output arguments for dense A with B>
819 %! [q, r, p] = qr (rand (3, 2), rand (3, 1));
820 %!error <too many output arguments for dense A with B>
821 %! [q, r, p] = qr (rand (3, 2), rand (3, 1), 0);
822 
823 %!function retval = __testqr (q, r, a, p)
824 %! tol = 100* eps (class (q));
825 %! retval = 0;
826 %! if (nargin == 3)
827 %! n1 = norm (q*r - a);
828 %! n2 = norm (q'*q - eye (columns (q)));
829 %! retval = (n1 < tol && n2 < tol);
830 %! else
831 %! n1 = norm (q'*q - eye (columns (q)));
832 %! retval = (n1 < tol);
833 %! if (isvector (p))
834 %! n2 = norm (q*r - a(:,p));
835 %! retval = (retval && n2 < tol);
836 %! else
837 %! n2 = norm (q*r - a*p);
838 %! retval = (retval && n2 < tol);
839 %! endif
840 %! endif
841 %!endfunction
842 
843 %!test
844 %! t = ones (24, 1);
845 %! j = 1;
846 %!
847 %! if (false) # eliminate big matrix tests
848 %! a = rand (5000, 20);
849 %! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
850 %! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
851 %! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
852 %! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
853 %!
854 %! a = a+1i*eps;
855 %! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
856 %! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
857 %! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
858 %! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
859 %! endif
860 %!
861 %! a = [ ones(1,15); sqrt(eps)*eye(15) ];
862 %! [q,r] = qr (a); t(j++) = __testqr (q, r, a);
863 %! [q,r] = qr (a'); t(j++) = __testqr (q, r, a');
864 %! [q,r,p] = qr (a); t(j++) = __testqr (q, r, a, p);
865 %! [q,r,p] = qr (a'); t(j++) = __testqr (q, r, a', p);
866 %!
867 %! a = a+1i*eps;
868 %! [q,r] = qr (a); t(j++) = __testqr (q, r, a);
869 %! [q,r] = qr (a'); t(j++) = __testqr (q, r, a');
870 %! [q,r,p] = qr (a); t(j++) = __testqr (q, r, a, p);
871 %! [q,r,p] = qr (a'); t(j++) = __testqr (q, r, a', p);
872 %!
873 %! a = [ ones(1,15); sqrt(eps)*eye(15) ];
874 %! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
875 %! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
876 %! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
877 %! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
878 %!
879 %! a = a+1i*eps;
880 %! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
881 %! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
882 %! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
883 %! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
884 %!
885 %! a = [ 611 196 -192 407 -8 -52 -49 29
886 %! 196 899 113 -192 -71 -43 -8 -44
887 %! -192 113 899 196 61 49 8 52
888 %! 407 -192 196 611 8 44 59 -23
889 %! -8 -71 61 8 411 -599 208 208
890 %! -52 -43 49 44 -599 411 208 208
891 %! -49 -8 8 59 208 208 99 -911
892 %! 29 -44 52 -23 208 208 -911 99 ];
893 %! [q,r] = qr (a);
894 %!
895 %! assert (all (t) && norm (q*r - a) < 5000*eps);
896 
897 %!test
898 %! a = single ([0, 2, 1; 2, 1, 2]);
899 %!
900 %! [q, r] = qr (a);
901 %! [qe, re] = qr (a, 0);
902 %!
903 %! assert (q * r, a, sqrt (eps ("single")));
904 %! assert (qe * re, a, sqrt (eps ("single")));
905 
906 %!test
907 %! a = single ([0, 2, 1; 2, 1, 2]);
908 %!
909 %! [q, r, p] = qr (a); # FIXME: not giving right dimensions.
910 %! [qe, re, pe] = qr (a, 0);
911 %!
912 %! assert (q * r, a * p, sqrt (eps ("single")));
913 %! assert (qe * re, a(:, pe), sqrt (eps ("single")));
914 
915 %!test
916 %! a = single ([0, 2; 2, 1; 1, 2]);
917 %!
918 %! [q, r] = qr (a);
919 %! [qe, re] = qr (a, 0);
920 %!
921 %! assert (q * r, a, sqrt (eps ("single")));
922 %! assert (qe * re, a, sqrt (eps ("single")));
923 
924 %!test
925 %! a = single ([0, 2; 2, 1; 1, 2]);
926 %!
927 %! [q, r, p] = qr (a);
928 %! [qe, re, pe] = qr (a, 0);
929 %!
930 %! assert (q * r, a * p, sqrt (eps ("single")));
931 %! assert (qe * re, a(:, pe), sqrt (eps ("single")));
932 
933 %!test
934 %! a = single ([0, 2, 1; 2, 1, 2; 3, 1, 2]);
935 %! b = single ([1, 3, 2; 1, 1, 0; 3, 0, 2]);
936 %!
937 %! [q, r] = qr (a);
938 %! [c, re] = qr (a, b);
939 %!
940 %! assert (r, re, sqrt (eps ("single")));
941 %! assert (q'*b, c, sqrt (eps ("single")));
942 
943 %!test
944 %! a = single ([0, 2, i; 2, 1, 2; 3, 1, 2]);
945 %! b = single ([1, 3, 2; 1, i, 0; 3, 0, 2]);
946 %!
947 %! [q, r] = qr (a);
948 %! [c, re] = qr (a, b);
949 %!
950 %! assert (r, re, sqrt (eps ("single")));
951 %! assert (q'*b, c, sqrt (eps ("single")));
952 
953 %!test
954 %! a = single ([0, 2, i; 2, 1, 2; 3, 1, 2]);
955 %! b = single ([1, 3, 2; 1, 1, 0; 3, 0, 2]);
956 %!
957 %! [q, r] = qr (a);
958 %! [c, re] = qr (a, b);
959 %!
960 %! assert (r, re, sqrt (eps));
961 %! assert (q'*b, c, sqrt (eps));
962 
963 %!test
964 %! a = single ([0, 2, 1; 2, 1, 2; 3, 1, 2]);
965 %! b = single ([1, 3, 2; 1, i, 0; 3, 0, 2]);
966 %!
967 %! [q, r] = qr (a);
968 %! [c, re] = qr (a, b);
969 %!
970 %! assert (r, re, sqrt (eps ("single")));
971 %! assert (q'*b, c, sqrt (eps ("single")));
972 
973 %!test
974 %! t = ones (24, 1);
975 %! j = 1;
976 %!
977 %! if (false) # eliminate big matrix tests
978 %! a = rand (5000,20);
979 %! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
980 %! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
981 %! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
982 %! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
983 %!
984 %! a = a+1i* eps ("single");
985 %! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
986 %! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
987 %! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
988 %! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
989 %! endif
990 %!
991 %! a = [ ones(1,15); sqrt(eps("single"))*eye(15) ];
992 %! [q,r] = qr (a); t(j++) = __testqr (q, r, a);
993 %! [q,r] = qr (a'); t(j++) = __testqr (q, r, a');
994 %! [q,r,p] = qr (a); t(j++) = __testqr (q, r, a, p);
995 %! [q,r,p] = qr (a'); t(j++) = __testqr (q, r, a', p);
996 %!
997 %! a = a+1i* eps ("single");
998 %! [q,r] = qr (a); t(j++) = __testqr (q, r, a);
999 %! [q,r] = qr (a'); t(j++) = __testqr (q, r, a');
1000 %! [q,r,p] = qr (a); t(j++) = __testqr (q, r, a, p);
1001 %! [q,r,p] = qr (a'); t(j++) = __testqr (q, r, a', p);
1002 %!
1003 %! a = [ ones(1,15); sqrt(eps("single"))*eye(15) ];
1004 %! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
1005 %! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
1006 %! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
1007 %! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
1008 %!
1009 %! a = a+1i* eps ("single");
1010 %! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
1011 %! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
1012 %! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
1013 %! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a',p);
1014 %!
1015 %! a = [ 611 196 -192 407 -8 -52 -49 29
1016 %! 196 899 113 -192 -71 -43 -8 -44
1017 %! -192 113 899 196 61 49 8 52
1018 %! 407 -192 196 611 8 44 59 -23
1019 %! -8 -71 61 8 411 -599 208 208
1020 %! -52 -43 49 44 -599 411 208 208
1021 %! -49 -8 8 59 208 208 99 -911
1022 %! 29 -44 52 -23 208 208 -911 99 ];
1023 %! [q,r] = qr (a);
1024 %!
1025 %! assert (all (t) && norm (q*r-a) < 5000* eps ("single"));
1026 
1027 ## The deactivated tests below can't be tested till rectangular back-subs is
1028 ## implemented for sparse matrices.
1029 
1030 %!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1031 %! n = 20; d = 0.2;
1032 %! ## initialize generators to make behavior reproducible
1033 %! rand ("state", 42);
1034 %! randn ("state", 42);
1035 %! a = sprandn (n,n,d) + speye (n,n);
1036 %! r = qr (a);
1037 %! assert (r'*r, a'*a, 1e-10);
1038 
1039 %!testif HAVE_COLAMD; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1040 %! n = 20; d = 0.2;
1041 %! ## initialize generators to make behavior reproducible
1042 %! rand ("state", 42);
1043 %! randn ("state", 42);
1044 %! a = sprandn (n,n,d) + speye (n,n);
1045 %! q = symamd (a);
1046 %! a = a(q,q);
1047 %! r = qr (a);
1048 %! assert (r'*r, a'*a, 1e-10);
1049 
1050 %!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1051 %! n = 20; d = 0.2;
1052 %! ## initialize generators to make behavior reproducible
1053 %! rand ("state", 42);
1054 %! randn ("state", 42);
1055 %! a = sprandn (n,n,d) + speye (n,n);
1056 %! [c,r] = qr (a, ones (n,1));
1057 %! assert (r\c, full (a)\ones (n,1), 10e-10);
1058 
1059 %!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1060 %! n = 20; d = 0.2;
1061 %! ## initialize generators to make behavior reproducible
1062 %! rand ("state", 42);
1063 %! randn ("state", 42);
1064 %! a = sprandn (n,n,d) + speye (n,n);
1065 %! b = randn (n,2);
1066 %! [c,r] = qr (a, b);
1067 %! assert (r\c, full (a)\b, 10e-10);
1068 
1069 ## Test under-determined systems!!
1070 %!#testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1071 %! n = 20; d = 0.2;
1072 %! ## initialize generators to make behavior reproducible
1073 %! rand ("state", 42);
1074 %! randn ("state", 42);
1075 %! a = sprandn (n,n+1,d) + speye (n,n+1);
1076 %! b = randn (n,2);
1077 %! [c,r] = qr (a, b);
1078 %! assert (r\c, full (a)\b, 10e-10);
1079 
1080 %!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1081 %! n = 20; d = 0.2;
1082 %! ## initialize generators to make behavior reproducible
1083 %! rand ("state", 42);
1084 %! randn ("state", 42);
1085 %! a = 1i* sprandn (n,n,d) + speye (n,n);
1086 %! r = qr (a);
1087 %! assert (r'*r,a'*a,1e-10);
1088 
1089 %!testif HAVE_COLAMD; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1090 %! n = 20; d = 0.2;
1091 %! ## initialize generators to make behavior reproducible
1092 %! rand ("state", 42);
1093 %! randn ("state", 42);
1094 %! a = 1i* sprandn (n,n,d) + speye (n,n);
1095 %! q = symamd (a);
1096 %! a = a(q,q);
1097 %! r = qr (a);
1098 %! assert (r'*r, a'*a, 1e-10);
1099 
1100 %!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1101 %! n = 20; d = 0.2;
1102 %! ## initialize generators to make behavior reproducible
1103 %! rand ("state", 42);
1104 %! randn ("state", 42);
1105 %! a = 1i* sprandn (n,n,d) + speye (n,n);
1106 %! [c,r] = qr (a, ones (n,1));
1107 %! assert (r\c, full (a)\ones (n,1), 10e-10);
1108 
1109 %!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1110 %! n = 20; d = 0.2;
1111 %! ## initialize generators to make behavior reproducible
1112 %! rand ("state", 42);
1113 %! randn ("state", 42);
1114 %! a = 1i* sprandn (n,n,d) + speye (n,n);
1115 %! b = randn (n,2);
1116 %! [c,r] = qr (a, b);
1117 %! assert (r\c, full (a)\b, 10e-10);
1118 
1119 ## Test under-determined systems!!
1120 %!#testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1121 %! n = 20; d = 0.2;
1122 %! ## initialize generators to make behavior reproducible
1123 %! rand ("state", 42);
1124 %! randn ("state", 42);
1125 %! a = 1i* sprandn (n,n+1,d) + speye (n,n+1);
1126 %! b = randn (n, 2);
1127 %! [c, r] = qr (a, b);
1128 %! assert (r\c, full (a)\b, 10e-10);
1129 
1130 %!testif HAVE_SPQR, HAVE_CHOLMOD
1131 %! n = 12; m = 20; d = 0.2;
1132 %! ## initialize generators to make behavior reproducible
1133 %! rand ("state", 42);
1134 %! randn ("state", 42);
1135 %! a = sprandn (m, n, d);
1136 %! b = randn (m, 2);
1137 %! [c, r] = qr (a, b);
1138 %! assert (r\c, full (a)\b, 10e-10);
1139 
1140 %!testif HAVE_SPQR, HAVE_CHOLMOD
1141 %! n = 12; m = 20; d = 0.2;
1142 %! ## initialize generators to make behavior reproducible
1143 %! rand ("state", 42);
1144 %! randn ("state", 42);
1145 %! a = sprandn (m, n, d);
1146 %! b = sprandn (m, 2, d);
1147 %! [c, r] = qr (a, b, 0);
1148 %! [c2, r2] = qr (full (a), full (b), 0);
1149 %! assert (r\c, r2\c2, 10e-10);
1150 
1151 %!testif HAVE_SPQR, HAVE_CHOLMOD
1152 %! n = 12; m = 20; d = 0.2;
1153 %! ## initialize generators to make behavior reproducible
1154 %! rand ("state", 42);
1155 %! randn ("state", 42);
1156 %! a = sprandn (m, n, d);
1157 %! b = randn (m, 2);
1158 %! [c, r, p] = qr (a, b, "matrix");
1159 %! x = p * (r\c);
1160 %! [c2, r2] = qr (full (a), b);
1161 %! x2 = r2 \ c2;
1162 %! assert (x, x2, 10e-10);
1163 
1164 %!testif HAVE_SPQR, HAVE_CHOLMOD
1165 %! n = 12; m = 20; d = 0.2;
1166 %! ## initialize generators to make behavior reproducible
1167 %! rand ("state", 42);
1168 %! randn ("state", 42);
1169 %! a = sprandn (m, n, d);
1170 %! [q, r, p] = qr (a, "matrix");
1171 %! assert (q * r, a * p, 10e-10);
1172 
1173 %!testif HAVE_SPQR, HAVE_CHOLMOD
1174 %! n = 12; m = 20; d = 0.2;
1175 %! ## initialize generators to make behavior reproducible
1176 %! rand ("state", 42);
1177 %! randn ("state", 42);
1178 %! a = sprandn (m, n, d);
1179 %! b = randn (m, 2);
1180 %! x = qr (a, b);
1181 %! [c2, r2] = qr (full (a), b);
1182 %! assert (x, r2\c2, 10e-10);
1183 
1184 %!testif HAVE_SPQR, HAVE_CHOLMOD
1185 %! n = 12; m = 20; d = 0.2;
1186 %! ## initialize generators to make behavior reproducible
1187 %! rand ("state", 42);
1188 %! randn ("state", 42);
1189 %! a = sprandn (m, n, d);
1190 %! b = i * randn (m, 2);
1191 %! x = qr (a, b);
1192 %! [c2, r2] = qr (full (a), b);
1193 %! assert (x, r2\c2, 10e-10);
1194 
1195 %!#testif HAVE_SPQR, HAVE_CHOLMOD
1196 %! n = 12; m = 20; d = 0.2;
1197 %! ## initialize generators to make behavior reproducible
1198 %! rand ("state", 42);
1199 %! randn ("state", 42);
1200 %! a = sprandn (m, n, d);
1201 %! b = i * randn (m, 2);
1202 %! [c, r] = qr (a, b);
1203 %! [c2, r2] = qr (full (a), b);
1204 %! assert (r\c, r2\c2, 10e-10);
1205 
1206 %!testif HAVE_SPQR, HAVE_CHOLMOD
1207 %! n = 12; m = 20; d = 0.2;
1208 %! ## initialize generators to make behavior reproducible
1209 %! rand ("state", 42);
1210 %! randn ("state", 42);
1211 %! a = sprandn (m, n, d);
1212 %! b = i * randn (m, 2);
1213 %! [c, r, p] = qr (a, b, "matrix");
1214 %! x = p * (r\c);
1215 %! [c2, r2] = qr (full (a), b);
1216 %! x2 = r2 \ c2;
1217 %! assert (x, x2, 10e-10);
1218 
1219 %!testif HAVE_SPQR, HAVE_CHOLMOD
1220 %! n = 12; m = 20; d = 0.2;
1221 %! ## initialize generators to make behavior reproducible
1222 %! rand ("state", 42);
1223 %! randn ("state", 42);
1224 %! a = i * sprandn (m, n, d);
1225 %! b = sprandn (m, 2, d);
1226 %! [c, r] = qr (a, b, 0);
1227 %! [c2, r2] = qr (full (a), full (b), 0);
1228 %! assert (r\c, r2\c2, 10e-10);
1229 
1230 %!testif HAVE_SPQR, HAVE_CHOLMOD
1231 %! n = 12; m = 20; d = 0.2;
1232 %! ## initialize generators to make behavior reproducible
1233 %! rand ("state", 42);
1234 %! randn ("state", 42);
1235 %! a = i * sprandn (m, n, d);
1236 %! b = randn (m, 2);
1237 %! [c, r, p] = qr (a, b, "matrix");
1238 %! x = p * (r\c);
1239 %! [c2, r2] = qr (full (a), b);
1240 %! x2 = r2 \ c2;
1241 %! assert(x, x2, 10e-10);
1242 
1243 %!testif HAVE_SPQR, HAVE_CHOLMOD
1244 %! n = 12; m = 20; d = 0.2;
1245 %! ## initialize generators to make behavior reproducible
1246 %! rand ("state", 42);
1247 %! randn ("state", 42);
1248 %! a = i * sprandn (m, n, d);
1249 %! [q, r, p] = qr (a, "matrix");
1250 %! assert(q * r, a * p, 10e-10);
1251 
1252 %!testif HAVE_SPQR, HAVE_CHOLMOD
1253 %! n = 12; m = 20; d = 0.2;
1254 %! ## initialize generators to make behavior reproducible
1255 %! rand ("state", 42);
1256 %! randn ("state", 42);
1257 %! a = i * sprandn (m, n, d);
1258 %! b = randn (m, 2);
1259 %! x = qr (a, b);
1260 %! [c2, r2] = qr (full (a), b);
1261 %! assert (x, r2\c2, 10e-10);
1262 
1263 %!testif HAVE_SPQR, HAVE_CHOLMOD
1264 %! a = sparse (5, 6);
1265 %! a(3,1) = 0.8;
1266 %! a(2,2) = 1.4;
1267 %! a(1,6) = -0.5;
1268 %! r = qr (a);
1269 %! assert (r'*r, a'*a, 10e-10);
1270 */
1271 
1272 static
1273 bool check_qr_dims (const octave_value& q, const octave_value& r,
1274  bool allow_ecf = false)
1275 {
1276  octave_idx_type m = q.rows ();
1277  octave_idx_type k = r.rows ();
1278  octave_idx_type n = r.columns ();
1279  return ((q.ndims () == 2 && r.ndims () == 2 && k == q.columns ())
1280  && (m == k || (allow_ecf && k == n && k < m)));
1281 }
1282 
1283 static
1284 bool check_index (const octave_value& i, bool vector_allowed = false)
1285 {
1286  return ((i.isreal () || i.isinteger ())
1287  && (i.is_scalar_type () || vector_allowed));
1288 }
1289 
1290 DEFUN (qrupdate, args, ,
1291  doc: /* -*- texinfo -*-
1292 @deftypefn {} {[@var{Q1}, @var{R1}] =} qrupdate (@var{Q}, @var{R}, @var{u}, @var{v})
1293 Update a QR factorization given update vectors or matrices.
1294 
1295 Given a QR@tie{}factorization of a real or complex matrix
1296 @w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and
1297 @var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization of
1298 @w{@var{A} + @var{u}*@var{v}'}, where @var{u} and @var{v} are column vectors
1299 (rank-1 update) or matrices with equal number of columns
1300 (rank-k update). Notice that the latter case is done as a sequence of
1301 rank-1 updates; thus, for k large enough, it will be both faster and more
1302 accurate to recompute the factorization from scratch.
1303 
1304 The QR@tie{}factorization supplied may be either full (Q is square) or
1305 economized (R is square).
1306 
1307 @seealso{qr, qrinsert, qrdelete, qrshift}
1308 @end deftypefn */)
1309 {
1310  octave_value_list retval;
1311 
1312  if (args.length () != 4)
1313  print_usage ();
1314 
1315  octave_value argq = args(0);
1316  octave_value argr = args(1);
1317  octave_value argu = args(2);
1318  octave_value argv = args(3);
1319 
1320  if (! argq.isnumeric () || ! argr.isnumeric ()
1321  || ! argu.isnumeric () || ! argv.isnumeric ())
1322  print_usage ();
1323 
1324  if (! check_qr_dims (argq, argr, true))
1325  error ("qrupdate: Q and R dimensions don't match");
1326 
1327  if (argq.isreal () && argr.isreal () && argu.isreal ()
1328  && argv.isreal ())
1329  {
1330  // all real case
1331  if (argq.is_single_type () || argr.is_single_type ()
1332  || argu.is_single_type () || argv.is_single_type ())
1333  {
1334  FloatMatrix Q = argq.float_matrix_value ();
1335  FloatMatrix R = argr.float_matrix_value ();
1336  FloatMatrix u = argu.float_matrix_value ();
1337  FloatMatrix v = argv.float_matrix_value ();
1338 
1339  math::qr<FloatMatrix> fact (Q, R);
1340  fact.update (u, v);
1341 
1342  retval = ovl (fact.Q (), get_qr_r (fact));
1343  }
1344  else
1345  {
1346  Matrix Q = argq.matrix_value ();
1347  Matrix R = argr.matrix_value ();
1348  Matrix u = argu.matrix_value ();
1349  Matrix v = argv.matrix_value ();
1350 
1351  math::qr<Matrix> fact (Q, R);
1352  fact.update (u, v);
1353 
1354  retval = ovl (fact.Q (), get_qr_r (fact));
1355  }
1356  }
1357  else
1358  {
1359  // complex case
1360  if (argq.is_single_type () || argr.is_single_type ()
1361  || argu.is_single_type () || argv.is_single_type ())
1362  {
1367 
1368  math::qr<FloatComplexMatrix> fact (Q, R);
1369  fact.update (u, v);
1370 
1371  retval = ovl (fact.Q (), get_qr_r (fact));
1372  }
1373  else
1374  {
1376  ComplexMatrix R = argr.complex_matrix_value ();
1377  ComplexMatrix u = argu.complex_matrix_value ();
1378  ComplexMatrix v = argv.complex_matrix_value ();
1379 
1380  math::qr<ComplexMatrix> fact (Q, R);
1381  fact.update (u, v);
1382 
1383  retval = ovl (fact.Q (), get_qr_r (fact));
1384  }
1385  }
1386 
1387  return retval;
1388 }
1389 
1390 /*
1391 %!shared A, u, v, Ac, uc, vc
1392 %! A = [0.091364 0.613038 0.999083;
1393 %! 0.594638 0.425302 0.603537;
1394 %! 0.383594 0.291238 0.085574;
1395 %! 0.265712 0.268003 0.238409;
1396 %! 0.669966 0.743851 0.445057 ];
1397 %!
1398 %! u = [0.85082;
1399 %! 0.76426;
1400 %! 0.42883;
1401 %! 0.53010;
1402 %! 0.80683 ];
1403 %!
1404 %! v = [0.98810;
1405 %! 0.24295;
1406 %! 0.43167 ];
1407 %!
1408 %! Ac = [0.620405 + 0.956953i 0.480013 + 0.048806i 0.402627 + 0.338171i;
1409 %! 0.589077 + 0.658457i 0.013205 + 0.279323i 0.229284 + 0.721929i;
1410 %! 0.092758 + 0.345687i 0.928679 + 0.241052i 0.764536 + 0.832406i;
1411 %! 0.912098 + 0.721024i 0.049018 + 0.269452i 0.730029 + 0.796517i;
1412 %! 0.112849 + 0.603871i 0.486352 + 0.142337i 0.355646 + 0.151496i ];
1413 %!
1414 %! uc = [0.20351 + 0.05401i;
1415 %! 0.13141 + 0.43708i;
1416 %! 0.29808 + 0.08789i;
1417 %! 0.69821 + 0.38844i;
1418 %! 0.74871 + 0.25821i ];
1419 %!
1420 %! vc = [0.85839 + 0.29468i;
1421 %! 0.20820 + 0.93090i;
1422 %! 0.86184 + 0.34689i ];
1423 %!
1424 
1425 %!test
1426 %! [Q,R] = qr (A);
1427 %! [Q,R] = qrupdate (Q, R, u, v);
1428 %! assert (norm (vec (Q'*Q - eye (5)), Inf) < 1e1*eps);
1429 %! assert (norm (vec (triu (R)-R), Inf) == 0);
1430 %! assert (norm (vec (Q*R - A - u*v'), Inf) < norm (A)*1e1*eps);
1431 %!
1432 %!test
1433 %! [Q,R] = qr (Ac);
1434 %! [Q,R] = qrupdate (Q, R, uc, vc);
1435 %! assert (norm (vec (Q'*Q - eye (5)), Inf) < 1e1*eps);
1436 %! assert (norm (vec (triu (R)-R), Inf) == 0);
1437 %! assert (norm (vec (Q*R - Ac - uc*vc'), Inf) < norm (Ac)*1e1*eps);
1438 
1439 %!test
1440 %! [Q,R] = qr (single (A));
1441 %! [Q,R] = qrupdate (Q, R, single (u), single (v));
1442 %! assert (norm (vec (Q'*Q - eye (5,"single")), Inf) < 1e1* eps ("single"));
1443 %! assert (norm (vec (triu (R)-R), Inf) == 0);
1444 %! assert (norm (vec (Q*R - single (A) - single (u)* single (v)'), Inf)
1445 %! < norm (single (A))*1e1* eps ("single"));
1446 %!
1447 %!test
1448 %! [Q,R] = qr (single (Ac));
1449 %! [Q,R] = qrupdate (Q, R, single (uc), single (vc));
1450 %! assert (norm (vec (Q'*Q - eye (5,"single")), Inf) < 1e1* eps ("single"));
1451 %! assert (norm (vec (triu (R)-R), Inf) == 0);
1452 %! assert (norm (vec (Q*R - single (Ac) - single (uc)* single (vc)'), Inf)
1453 %! < norm (single (Ac))*1e1* eps ("single"));
1454 */
1455 
1456 DEFUN (qrinsert, args, ,
1457  doc: /* -*- texinfo -*-
1458 @deftypefn {} {[@var{Q1}, @var{R1}] =} qrinsert (@var{Q}, @var{R}, @var{j}, @var{x}, @var{orient})
1459 Update a QR factorization given a row or column to insert in the original
1460 factored matrix.
1461 
1462 
1463 Given a QR@tie{}factorization of a real or complex matrix
1464 @w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and
1465 @var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization of
1466 @w{[A(:,1:j-1) x A(:,j:n)]}, where @var{u} is a column vector to be inserted
1467 into @var{A} (if @var{orient} is @qcode{"col"}), or the
1468 QR@tie{}factorization of @w{[A(1:j-1,:);x;A(:,j:n)]}, where @var{x} is a row
1469 vector to be inserted into @var{A} (if @var{orient} is @qcode{"row"}).
1470 
1471 The default value of @var{orient} is @qcode{"col"}. If @var{orient} is
1472 @qcode{"col"}, @var{u} may be a matrix and @var{j} an index vector
1473 resulting in the QR@tie{}factorization of a matrix @var{B} such that
1474 @w{B(:,@var{j})} gives @var{u} and @w{B(:,@var{j}) = []} gives @var{A}.
1475 Notice that the latter case is done as a sequence of k insertions;
1476 thus, for k large enough, it will be both faster and more accurate to
1477 recompute the factorization from scratch.
1478 
1479 If @var{orient} is @qcode{"col"}, the QR@tie{}factorization supplied may
1480 be either full (Q is square) or economized (R is square).
1481 
1482 If @var{orient} is @qcode{"row"}, full factorization is needed.
1483 @seealso{qr, qrupdate, qrdelete, qrshift}
1484 @end deftypefn */)
1485 {
1486  int nargin = args.length ();
1487 
1488  if (nargin < 4 || nargin > 5)
1489  print_usage ();
1490 
1491  octave_value argq = args(0);
1492  octave_value argr = args(1);
1493  octave_value argj = args(2);
1494  octave_value argx = args(3);
1495 
1496  if (! argq.isnumeric () || ! argr.isnumeric ()
1497  || ! argx.isnumeric ()
1498  || (nargin > 4 && ! args(4).is_string ()))
1499  print_usage ();
1500 
1501  std::string orient = (nargin < 5) ? "col" : args(4).string_value ();
1502  bool col = (orient == "col");
1503 
1504  if (! col && orient != "row")
1505  error (R"(qrinsert: ORIENT must be "col" or "row")");
1506 
1507  if (! check_qr_dims (argq, argr, col) || (! col && argx.rows () != 1))
1508  error ("qrinsert: dimension mismatch");
1509 
1510  if (! check_index (argj, col))
1511  error ("qrinsert: invalid index J");
1512 
1513  octave_value_list retval;
1514 
1516 
1517  octave_idx_type one = 1;
1518 
1519  if (argq.isreal () && argr.isreal () && argx.isreal ())
1520  {
1521  // real case
1522  if (argq.is_single_type () || argr.is_single_type ()
1523  || argx.is_single_type ())
1524  {
1525  FloatMatrix Q = argq.float_matrix_value ();
1526  FloatMatrix R = argr.float_matrix_value ();
1527  FloatMatrix x = argx.float_matrix_value ();
1528 
1529  math::qr<FloatMatrix> fact (Q, R);
1530 
1531  if (col)
1532  fact.insert_col (x, j-one);
1533  else
1534  fact.insert_row (x.row (0), j(0)-one);
1535 
1536  retval = ovl (fact.Q (), get_qr_r (fact));
1537  }
1538  else
1539  {
1540  Matrix Q = argq.matrix_value ();
1541  Matrix R = argr.matrix_value ();
1542  Matrix x = argx.matrix_value ();
1543 
1544  math::qr<Matrix> fact (Q, R);
1545 
1546  if (col)
1547  fact.insert_col (x, j-one);
1548  else
1549  fact.insert_row (x.row (0), j(0)-one);
1550 
1551  retval = ovl (fact.Q (), get_qr_r (fact));
1552  }
1553  }
1554  else
1555  {
1556  // complex case
1557  if (argq.is_single_type () || argr.is_single_type ()
1558  || argx.is_single_type ())
1559  {
1563 
1564  math::qr<FloatComplexMatrix> fact (Q, R);
1565 
1566  if (col)
1567  fact.insert_col (x, j-one);
1568  else
1569  fact.insert_row (x.row (0), j(0)-one);
1570 
1571  retval = ovl (fact.Q (), get_qr_r (fact));
1572  }
1573  else
1574  {
1576  ComplexMatrix R = argr.complex_matrix_value ();
1578 
1579  math::qr<ComplexMatrix> fact (Q, R);
1580 
1581  if (col)
1582  fact.insert_col (x, j-one);
1583  else
1584  fact.insert_row (x.row (0), j(0)-one);
1585 
1586  retval = ovl (fact.Q (), get_qr_r (fact));
1587  }
1588  }
1589 
1590  return retval;
1591 }
1592 
1593 /*
1594 %!test
1595 %! [Q,R] = qr (A);
1596 %! [Q,R] = qrinsert (Q, R, 3, u);
1597 %! assert (norm (vec (Q'*Q - eye (5)), Inf) < 1e1*eps);
1598 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1599 %! assert (norm (vec (Q*R - [A(:,1:2) u A(:,3)]), Inf) < norm (A)*1e1*eps);
1600 %!test
1601 %! [Q,R] = qr (Ac);
1602 %! [Q,R] = qrinsert (Q, R, 3, uc);
1603 %! assert (norm (vec (Q'*Q - eye (5)), Inf) < 1e1*eps);
1604 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1605 %! assert (norm (vec (Q*R - [Ac(:,1:2) uc Ac(:,3)]), Inf) < norm (Ac)*1e1*eps);
1606 %!test
1607 %! x = [0.85082 0.76426 0.42883 ];
1608 %!
1609 %! [Q,R] = qr (A);
1610 %! [Q,R] = qrinsert (Q, R, 3, x, "row");
1611 %! assert (norm (vec (Q'*Q - eye (6)), Inf) < 1e1*eps);
1612 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1613 %! assert (norm (vec (Q*R - [A(1:2,:);x;A(3:5,:)]), Inf) < norm (A)*1e1*eps);
1614 %!test
1615 %! x = [0.20351 + 0.05401i 0.13141 + 0.43708i 0.29808 + 0.08789i ];
1616 %!
1617 %! [Q,R] = qr (Ac);
1618 %! [Q,R] = qrinsert (Q, R, 3, x, "row");
1619 %! assert (norm (vec (Q'*Q - eye (6)), Inf) < 1e1*eps);
1620 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1621 %! assert (norm (vec (Q*R - [Ac(1:2,:);x;Ac(3:5,:)]), Inf) < norm (Ac)*1e1*eps);
1622 
1623 %!test
1624 %! [Q,R] = qr (single (A));
1625 %! [Q,R] = qrinsert (Q, R, 3, single (u));
1626 %! assert (norm (vec (Q'*Q - eye (5,"single")), Inf) < 1e1* eps ("single"));
1627 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1628 %! assert (norm (vec (Q*R - single ([A(:,1:2) u A(:,3)])), Inf)
1629 %! < norm (single (A))*1e1* eps ("single"));
1630 %!test
1631 %! [Q,R] = qr (single (Ac));
1632 %! [Q,R] = qrinsert (Q, R, 3, single (uc));
1633 %! assert (norm (vec (Q'*Q - eye (5,"single")), Inf) < 1e1* eps ("single"));
1634 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1635 %! assert (norm (vec (Q*R - single ([Ac(:,1:2) uc Ac(:,3)])), Inf)
1636 %! < norm (single (Ac))*1e1* eps ("single"));
1637 %!test
1638 %! x = single ([0.85082 0.76426 0.42883 ]);
1639 %!
1640 %! [Q,R] = qr (single (A));
1641 %! [Q,R] = qrinsert (Q, R, 3, x, "row");
1642 %! assert (norm (vec (Q'*Q - eye (6,"single")), Inf) < 1e1* eps ("single"));
1643 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1644 %! assert (norm (vec (Q*R - single ([A(1:2,:);x;A(3:5,:)])), Inf)
1645 %! < norm (single (A))*1e1* eps ("single"));
1646 %!test
1647 %! x = single ([0.20351 + 0.05401i 0.13141 + 0.43708i 0.29808 + 0.08789i ]);
1648 %!
1649 %! [Q,R] = qr (single (Ac));
1650 %! [Q,R] = qrinsert (Q, R, 3, x, "row");
1651 %! assert (norm (vec (Q'*Q - eye (6,"single")), Inf) < 1e1* eps ("single"));
1652 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1653 %! assert (norm (vec (Q*R - single ([Ac(1:2,:);x;Ac(3:5,:)])), Inf)
1654 %! < norm (single (Ac))*1e1* eps ("single"));
1655 */
1656 
1657 DEFUN (qrdelete, args, ,
1658  doc: /* -*- texinfo -*-
1659 @deftypefn {} {[@var{Q1}, @var{R1}] =} qrdelete (@var{Q}, @var{R}, @var{j}, @var{orient})
1660 Update a QR factorization given a row or column to delete from the original
1661 factored matrix.
1662 
1663 Given a QR@tie{}factorization of a real or complex matrix
1664 @w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and
1665 @var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization of
1666 @w{[A(:,1:j-1), U, A(:,j:n)]},
1667 where @var{u} is a column vector to be inserted into @var{A}
1668 (if @var{orient} is @qcode{"col"}),
1669 or the QR@tie{}factorization of @w{[A(1:j-1,:);X;A(:,j:n)]},
1670 where @var{x} is a row @var{orient} is @qcode{"row"}).
1671 The default value of @var{orient} is @qcode{"col"}.
1672 
1673 If @var{orient} is @qcode{"col"}, @var{j} may be an index vector
1674 resulting in the QR@tie{}factorization of a matrix @var{B} such that
1675 @w{A(:,@var{j}) = []} gives @var{B}. Notice that the latter case is done as
1676 a sequence of k deletions; thus, for k large enough, it will be both faster
1677 and more accurate to recompute the factorization from scratch.
1678 
1679 If @var{orient} is @qcode{"col"}, the QR@tie{}factorization supplied may
1680 be either full (Q is square) or economized (R is square).
1681 
1682 If @var{orient} is @qcode{"row"}, full factorization is needed.
1683 @seealso{qr, qrupdate, qrinsert, qrshift}
1684 @end deftypefn */)
1685 {
1686  int nargin = args.length ();
1687 
1688  if (nargin < 3 || nargin > 4)
1689  print_usage ();
1690 
1691  octave_value argq = args(0);
1692  octave_value argr = args(1);
1693  octave_value argj = args(2);
1694 
1695  if (! argq.isnumeric () || ! argr.isnumeric ()
1696  || (nargin > 3 && ! args(3).is_string ()))
1697  print_usage ();
1698 
1699  std::string orient = (nargin < 4) ? "col" : args(3).string_value ();
1700  bool col = orient == "col";
1701 
1702  if (! col && orient != "row")
1703  error (R"(qrdelete: ORIENT must be "col" or "row")");
1704 
1705  if (! check_qr_dims (argq, argr, col))
1706  error ("qrdelete: dimension mismatch");
1707 
1709  if (! check_index (argj, col))
1710  error ("qrdelete: invalid index J");
1711 
1712  octave_value_list retval;
1713 
1714  octave_idx_type one = 1;
1715 
1716  if (argq.isreal () && argr.isreal ())
1717  {
1718  // real case
1719  if (argq.is_single_type () || argr.is_single_type ())
1720  {
1721  FloatMatrix Q = argq.float_matrix_value ();
1722  FloatMatrix R = argr.float_matrix_value ();
1723 
1724  math::qr<FloatMatrix> fact (Q, R);
1725 
1726  if (col)
1727  fact.delete_col (j-one);
1728  else
1729  fact.delete_row (j(0)-one);
1730 
1731  retval = ovl (fact.Q (), get_qr_r (fact));
1732  }
1733  else
1734  {
1735  Matrix Q = argq.matrix_value ();
1736  Matrix R = argr.matrix_value ();
1737 
1738  math::qr<Matrix> fact (Q, R);
1739 
1740  if (col)
1741  fact.delete_col (j-one);
1742  else
1743  fact.delete_row (j(0)-one);
1744 
1745  retval = ovl (fact.Q (), get_qr_r (fact));
1746  }
1747  }
1748  else
1749  {
1750  // complex case
1751  if (argq.is_single_type () || argr.is_single_type ())
1752  {
1755 
1756  math::qr<FloatComplexMatrix> fact (Q, R);
1757 
1758  if (col)
1759  fact.delete_col (j-one);
1760  else
1761  fact.delete_row (j(0)-one);
1762 
1763  retval = ovl (fact.Q (), get_qr_r (fact));
1764  }
1765  else
1766  {
1768  ComplexMatrix R = argr.complex_matrix_value ();
1769 
1770  math::qr<ComplexMatrix> fact (Q, R);
1771 
1772  if (col)
1773  fact.delete_col (j-one);
1774  else
1775  fact.delete_row (j(0)-one);
1776 
1777  retval = ovl (fact.Q (), get_qr_r (fact));
1778  }
1779  }
1780 
1781  return retval;
1782 }
1783 
1784 /*
1785 %!test
1786 %! AA = [0.091364 0.613038 0.027504 0.999083;
1787 %! 0.594638 0.425302 0.562834 0.603537;
1788 %! 0.383594 0.291238 0.742073 0.085574;
1789 %! 0.265712 0.268003 0.783553 0.238409;
1790 %! 0.669966 0.743851 0.457255 0.445057 ];
1791 %!
1792 %! [Q,R] = qr (AA);
1793 %! [Q,R] = qrdelete (Q, R, 3);
1794 %! assert (norm (vec (Q'*Q - eye (5)), Inf) < 16*eps);
1795 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1796 %! assert (norm (vec (Q*R - [AA(:,1:2) AA(:,4)]), Inf) < norm (AA)*1e1*eps);
1797 %!
1798 %!test
1799 %! AA = [0.364554 + 0.993117i 0.669818 + 0.510234i 0.426568 + 0.041337i 0.847051 + 0.233291i;
1800 %! 0.049600 + 0.242783i 0.448946 + 0.484022i 0.141155 + 0.074420i 0.446746 + 0.392706i;
1801 %! 0.581922 + 0.657416i 0.581460 + 0.030016i 0.219909 + 0.447288i 0.201144 + 0.069132i;
1802 %! 0.694986 + 0.000571i 0.682327 + 0.841712i 0.807537 + 0.166086i 0.192767 + 0.358098i;
1803 %! 0.945002 + 0.066788i 0.350492 + 0.642638i 0.579629 + 0.048102i 0.600170 + 0.636938i ] * I;
1804 %!
1805 %! [Q,R] = qr (AA);
1806 %! [Q,R] = qrdelete (Q, R, 3);
1807 %! assert (norm (vec (Q'*Q - eye (5)), Inf) < 16*eps);
1808 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1809 %! assert (norm (vec (Q*R - [AA(:,1:2) AA(:,4)]), Inf) < norm (AA)*1e1*eps);
1810 %!
1811 %!test
1812 %! AA = [0.091364 0.613038 0.027504 0.999083;
1813 %! 0.594638 0.425302 0.562834 0.603537;
1814 %! 0.383594 0.291238 0.742073 0.085574;
1815 %! 0.265712 0.268003 0.783553 0.238409;
1816 %! 0.669966 0.743851 0.457255 0.445057 ];
1817 %!
1818 %! [Q,R] = qr (AA);
1819 %! [Q,R] = qrdelete (Q, R, 3, "row");
1820 %! assert (norm (vec (Q'*Q - eye (4)), Inf) < 1e1*eps);
1821 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1822 %! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf) < norm (AA)*1e1*eps);
1823 %!
1824 %!test
1825 %! AA = [0.364554 + 0.993117i 0.669818 + 0.510234i 0.426568 + 0.041337i 0.847051 + 0.233291i;
1826 %! 0.049600 + 0.242783i 0.448946 + 0.484022i 0.141155 + 0.074420i 0.446746 + 0.392706i;
1827 %! 0.581922 + 0.657416i 0.581460 + 0.030016i 0.219909 + 0.447288i 0.201144 + 0.069132i;
1828 %! 0.694986 + 0.000571i 0.682327 + 0.841712i 0.807537 + 0.166086i 0.192767 + 0.358098i;
1829 %! 0.945002 + 0.066788i 0.350492 + 0.642638i 0.579629 + 0.048102i 0.600170 + 0.636938i ] * I;
1830 %!
1831 %! [Q,R] = qr (AA);
1832 %! [Q,R] = qrdelete (Q, R, 3, "row");
1833 %! assert (norm (vec (Q'*Q - eye (4)), Inf) < 1e1*eps);
1834 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1835 %! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf) < norm (AA)*1e1*eps);
1836 
1837 %!test
1838 %! AA = single ([0.091364 0.613038 0.027504 0.999083;
1839 %! 0.594638 0.425302 0.562834 0.603537;
1840 %! 0.383594 0.291238 0.742073 0.085574;
1841 %! 0.265712 0.268003 0.783553 0.238409;
1842 %! 0.669966 0.743851 0.457255 0.445057 ]);
1843 %!
1844 %! [Q,R] = qr (AA);
1845 %! [Q,R] = qrdelete (Q, R, 3);
1846 %! assert (norm (vec (Q'*Q - eye (5,"single")), Inf) < 1e1* eps ("single"));
1847 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1848 %! assert (norm (vec (Q*R - [AA(:,1:2) AA(:,4)]), Inf)
1849 %! < norm (AA)*1e1* eps ("single"));
1850 %!
1851 %!test
1852 %! AA = single ([0.364554 + 0.993117i 0.669818 + 0.510234i 0.426568 + 0.041337i 0.847051 + 0.233291i;
1853 %! 0.049600 + 0.242783i 0.448946 + 0.484022i 0.141155 + 0.074420i 0.446746 + 0.392706i;
1854 %! 0.581922 + 0.657416i 0.581460 + 0.030016i 0.219909 + 0.447288i 0.201144 + 0.069132i;
1855 %! 0.694986 + 0.000571i 0.682327 + 0.841712i 0.807537 + 0.166086i 0.192767 + 0.358098i;
1856 %! 0.945002 + 0.066788i 0.350492 + 0.642638i 0.579629 + 0.048102i 0.600170 + 0.636938i ]) * I;
1857 %!
1858 %! [Q,R] = qr (AA);
1859 %! [Q,R] = qrdelete (Q, R, 3);
1860 %! assert (norm (vec (Q'*Q - eye (5,"single")), Inf) < 1e1* eps ("single"));
1861 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1862 %! assert (norm (vec (Q*R - [AA(:,1:2) AA(:,4)]), Inf)
1863 %! < norm (AA)*1e1* eps ("single"));
1864 
1865 %!test
1866 %! AA = single ([0.091364 0.613038 0.027504 0.999083;
1867 %! 0.594638 0.425302 0.562834 0.603537;
1868 %! 0.383594 0.291238 0.742073 0.085574;
1869 %! 0.265712 0.268003 0.783553 0.238409;
1870 %! 0.669966 0.743851 0.457255 0.445057 ]);
1871 %!
1872 %! [Q,R] = qr (AA);
1873 %! [Q,R] = qrdelete (Q, R, 3, "row");
1874 %! assert (norm (vec (Q'*Q - eye (4,"single")), Inf) < 1.5e1* eps ("single"));
1875 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1876 %! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf)
1877 %! < norm (AA)*1e1* eps ("single"));
1878 %!testif HAVE_QRUPDATE
1879 %! ## Same test as above but with more precicision
1880 %! AA = single ([0.091364 0.613038 0.027504 0.999083;
1881 %! 0.594638 0.425302 0.562834 0.603537;
1882 %! 0.383594 0.291238 0.742073 0.085574;
1883 %! 0.265712 0.268003 0.783553 0.238409;
1884 %! 0.669966 0.743851 0.457255 0.445057 ]);
1885 %!
1886 %! [Q,R] = qr (AA);
1887 %! [Q,R] = qrdelete (Q, R, 3, "row");
1888 %! assert (norm (vec (Q'*Q - eye (4,"single")), Inf) < 1e1* eps ("single"));
1889 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1890 %! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf)
1891 %! < norm (AA)*1e1* eps ("single"));
1892 %!
1893 %!test
1894 %! AA = single ([0.364554 + 0.993117i 0.669818 + 0.510234i 0.426568 + 0.041337i 0.847051 + 0.233291i;
1895 %! 0.049600 + 0.242783i 0.448946 + 0.484022i 0.141155 + 0.074420i 0.446746 + 0.392706i;
1896 %! 0.581922 + 0.657416i 0.581460 + 0.030016i 0.219909 + 0.447288i 0.201144 + 0.069132i;
1897 %! 0.694986 + 0.000571i 0.682327 + 0.841712i 0.807537 + 0.166086i 0.192767 + 0.358098i;
1898 %! 0.945002 + 0.066788i 0.350492 + 0.642638i 0.579629 + 0.048102i 0.600170 + 0.636938i ]) * I;
1899 %!
1900 %! [Q,R] = qr (AA);
1901 %! [Q,R] = qrdelete (Q, R, 3, "row");
1902 %! assert (norm (vec (Q'*Q - eye (4,"single")), Inf) < 1e1* eps ("single"));
1903 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1904 %! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf)
1905 %! < norm (AA)*1e1* eps ("single"));
1906 */
1907 
1908 DEFUN (qrshift, args, ,
1909  doc: /* -*- texinfo -*-
1910 @deftypefn {} {[@var{Q1}, @var{R1}] =} qrshift (@var{Q}, @var{R}, @var{i}, @var{j})
1911 Update a QR factorization given a range of columns to shift in the original
1912 factored matrix.
1913 
1914 Given a QR@tie{}factorization of a real or complex matrix
1915 @w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and
1916 @var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization
1917 of @w{@var{A}(:,p)}, where @w{p} is the permutation @*
1918 @code{p = [1:i-1, shift(i:j, 1), j+1:n]} if @w{@var{i} < @var{j}} @*
1919  or @*
1920 @code{p = [1:j-1, shift(j:i,-1), i+1:n]} if @w{@var{j} < @var{i}}. @*
1921 
1922 @seealso{qr, qrupdate, qrinsert, qrdelete}
1923 @end deftypefn */)
1924 {
1925  if (args.length () != 4)
1926  print_usage ();
1927 
1928  octave_value argq = args(0);
1929  octave_value argr = args(1);
1930  octave_value argi = args(2);
1931  octave_value argj = args(3);
1932 
1933  if (! argq.isnumeric () || ! argr.isnumeric ())
1934  print_usage ();
1935 
1936  if (! check_qr_dims (argq, argr, true))
1937  error ("qrshift: dimensions mismatch");
1938 
1939  octave_idx_type i = argi.idx_type_value ();
1940  octave_idx_type j = argj.idx_type_value ();
1941 
1942  if (! check_index (argi) || ! check_index (argj))
1943  error ("qrshift: invalid index I or J");
1944 
1945  octave_value_list retval;
1946 
1947  if (argq.isreal () && argr.isreal ())
1948  {
1949  // all real case
1950  if (argq.is_single_type ()
1951  && argr.is_single_type ())
1952  {
1953  FloatMatrix Q = argq.float_matrix_value ();
1954  FloatMatrix R = argr.float_matrix_value ();
1955 
1956  math::qr<FloatMatrix> fact (Q, R);
1957  fact.shift_cols (i-1, j-1);
1958 
1959  retval = ovl (fact.Q (), get_qr_r (fact));
1960  }
1961  else
1962  {
1963  Matrix Q = argq.matrix_value ();
1964  Matrix R = argr.matrix_value ();
1965 
1966  math::qr<Matrix> fact (Q, R);
1967  fact.shift_cols (i-1, j-1);
1968 
1969  retval = ovl (fact.Q (), get_qr_r (fact));
1970  }
1971  }
1972  else
1973  {
1974  // complex case
1975  if (argq.is_single_type ()
1976  && argr.is_single_type ())
1977  {
1980 
1981  math::qr<FloatComplexMatrix> fact (Q, R);
1982  fact.shift_cols (i-1, j-1);
1983 
1984  retval = ovl (fact.Q (), get_qr_r (fact));
1985  }
1986  else
1987  {
1989  ComplexMatrix R = argr.complex_matrix_value ();
1990 
1991  math::qr<ComplexMatrix> fact (Q, R);
1992  fact.shift_cols (i-1, j-1);
1993 
1994  retval = ovl (fact.Q (), get_qr_r (fact));
1995  }
1996  }
1997 
1998  return retval;
1999 }
2000 
2001 /*
2002 %!test
2003 %! AA = A.';
2004 %! i = 2; j = 4; p = [1:i-1, shift(i:j,-1), j+1:5];
2005 %!
2006 %! [Q,R] = qr (AA);
2007 %! [Q,R] = qrshift (Q, R, i, j);
2008 %! assert (norm (vec (Q'*Q - eye (3)), Inf) < 1e1*eps);
2009 %! assert (norm (vec (triu (R) - R), Inf) == 0);
2010 %! assert (norm (vec (Q*R - AA(:,p)), Inf) < norm (AA)*1e1*eps);
2011 %!
2012 %! j = 2; i = 4; p = [1:j-1, shift(j:i,+1), i+1:5];
2013 %!
2014 %! [Q,R] = qr (AA);
2015 %! [Q,R] = qrshift (Q, R, i, j);
2016 %! assert (norm (vec (Q'*Q - eye (3)), Inf) < 1e1*eps);
2017 %! assert (norm (vec (triu (R) - R), Inf) == 0);
2018 %! assert (norm (vec (Q*R - AA(:,p)), Inf) < norm (AA)*1e1*eps);
2019 %!
2020 %!test
2021 %! AA = Ac.';
2022 %! i = 2; j = 4; p = [1:i-1, shift(i:j,-1), j+1:5];
2023 %!
2024 %! [Q,R] = qr (AA);
2025 %! [Q,R] = qrshift (Q, R, i, j);
2026 %! assert (norm (vec (Q'*Q - eye (3)), Inf) < 1e1*eps);
2027 %! assert (norm (vec (triu (R) - R), Inf) == 0);
2028 %! assert (norm (vec (Q*R - AA(:,p)), Inf) < norm (AA)*1e1*eps);
2029 %!
2030 %! j = 2; i = 4; p = [1:j-1, shift(j:i,+1), i+1:5];
2031 %!
2032 %! [Q,R] = qr (AA);
2033 %! [Q,R] = qrshift (Q, R, i, j);
2034 %! assert (norm (vec (Q'*Q - eye (3)), Inf) < 1e1*eps);
2035 %! assert (norm (vec (triu (R) - R), Inf) == 0);
2036 %! assert (norm (vec (Q*R - AA(:,p)), Inf) < norm (AA)*1e1*eps);
2037 
2038 %!test
2039 %! AA = single (A).';
2040 %! i = 2; j = 4; p = [1:i-1, shift(i:j,-1), j+1:5];
2041 %!
2042 %! [Q,R] = qr (AA);
2043 %! [Q,R] = qrshift (Q, R, i, j);
2044 %! assert (norm (vec (Q'*Q - eye (3,"single")), Inf) < 1e1* eps ("single"));
2045 %! assert (norm (vec (triu (R) - R), Inf) == 0);
2046 %! assert (norm (vec (Q*R - AA(:,p)), Inf) < norm (AA)*1e1* eps ("single"));
2047 %!
2048 %! j = 2; i = 4; p = [1:j-1, shift(j:i,+1), i+1:5];
2049 %!
2050 %! [Q,R] = qr (AA);
2051 %! [Q,R] = qrshift (Q, R, i, j);
2052 %! assert (norm (vec (Q'*Q - eye (3,"single")), Inf) < 1e1* eps ("single"));
2053 %! assert (norm (vec (triu (R) - R), Inf) == 0);
2054 %! assert (norm (vec (Q*R - AA(:,p)), Inf) < norm (AA)*1e1* eps ("single"));
2055 %!
2056 %!test
2057 %! AA = single (Ac).';
2058 %! i = 2; j = 4; p = [1:i-1, shift(i:j,-1), j+1:5];
2059 %!
2060 %! [Q,R] = qr (AA);
2061 %! [Q,R] = qrshift (Q, R, i, j);
2062 %! assert (norm (vec (Q'*Q - eye (3,"single")), Inf) < 1e1* eps ("single"));
2063 %! assert (norm (vec (triu (R) - R), Inf) == 0);
2064 %! assert (norm (vec (Q*R - AA(:,p)), Inf) < norm (AA)*1e1* eps ("single"));
2065 %!
2066 %! j = 2; i = 4; p = [1:j-1, shift(j:i,+1), i+1:5];
2067 %!
2068 %! [Q,R] = qr (AA);
2069 %! [Q,R] = qrshift (Q, R, i, j);
2070 %! assert (norm (vec (Q'*Q - eye (3,"single")), Inf) < 1e1* eps ("single"));
2071 %! assert (norm (vec (triu (R) - R), Inf) == 0);
2072 %! assert (norm (vec (Q*R - AA(:,p)), Inf) < norm (AA)*1e1* eps ("single"));
2073 */
2074 
OCTAVE_END_NAMESPACE(octave)
ComplexColumnVector conj(const ComplexColumnVector &a)
Definition: CColVector.cc:217
Definition: dMatrix.h:42
SparseMatrix sparse_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:945
bool isreal(void) const
Definition: ov.h:783
bool issparse(void) const
Definition: ov.h:798
ComplexMatrix complex_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:916
octave_idx_type rows(void) const
Definition: ov.h:590
bool isnumeric(void) const
Definition: ov.h:795
octave_idx_type idx_type_value(bool req_int=false, bool frc_str_conv=false) const
bool is_scalar_type(void) const
Definition: ov.h:789
bool isinteger(void) const
Definition: ov.h:775
octave_idx_type columns(void) const
Definition: ov.h:592
int ndims(void) const
Definition: ov.h:596
FloatMatrix float_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:901
bool is_single_type(void) const
Definition: ov.h:743
Matrix matrix_value(bool frc_str_conv=false) const
Definition: ov.h:898
FloatComplexMatrix float_complex_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:920
bool iscomplex(void) const
Definition: ov.h:786
OCTINTERP_API Array< octave_idx_type > octave_idx_type_vector_value(bool req_int=false, bool frc_str_conv=false, bool frc_vec_conv=false) const
SparseComplexMatrix sparse_complex_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:949
Definition: qr.h:40
OCTAVE_BEGIN_NAMESPACE(octave) static octave_value daspk_fcn
OCTINTERP_API void print_usage(void)
Definition: defun-int.h:72
#define DEFUN(name, args_name, nargout_name, doc)
Macro to define a builtin function.
Definition: defun.h:56
void error(const char *fmt,...)
Definition: error.cc:979
void err_wrong_type_arg(const char *name, const char *s)
Definition: errwarn.cc:166
static bool check_index(const octave_value &i, bool vector_allowed=false)
Definition: qr.cc:1284
static bool check_qr_dims(const octave_value &q, const octave_value &r, bool allow_ecf=false)
Definition: qr.cc:1273
static math::qr< T >::type qr_type(int nargout, bool economy)
Definition: qr.cc:74
static octave_value get_qr_r(const math::qr< MT > &fact)
Definition: qr.cc:63
F77_RET_T const F77_INT const F77_INT const F77_INT F77_DBLE const F77_INT F77_DBLE const F77_INT F77_DBLE * Q
F77_RET_T const F77_DBLE * x
T octave_idx_type m
Definition: mx-inlines.cc:773
octave_idx_type n
Definition: mx-inlines.cc:753
T * r
Definition: mx-inlines.cc:773
return octave_value(v1.char_array_value() . concat(v2.char_array_value(), ra_idx),((a1.is_sq_string()||a2.is_sq_string()) ? '\'' :'"'))
octave_value_list ovl(const OV_Args &... args)
Construct an octave_value_list with less typing.
Definition: ovl.h:211
template class OCTAVE_EXTERN_TEMPLATE_API sparse_qr< SparseMatrix >
template class OCTAVE_EXTERN_TEMPLATE_API sparse_qr< SparseComplexMatrix >