GNU Octave  9.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
qr.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 //
3 // Copyright (C) 1996-2024 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  {
407  math::sparse_qr<SparseComplexMatrix>
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  {
414  math::sparse_qr<SparseComplexMatrix>
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  {
427  math::sparse_qr<SparseComplexMatrix>
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  {
437  math::sparse_qr<SparseComplexMatrix>
438  q (arg.sparse_complex_matrix_value (), 0);
439  retval = ovl (q.Q (economy), q.R (economy));
440  }
441  else
442  {
443  math::sparse_qr<SparseComplexMatrix>
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  {
480  math::sparse_qr<SparseMatrix>
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  {
487  math::sparse_qr<SparseMatrix>
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  {
501  math::sparse_qr<SparseMatrix>
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  {
511  math::sparse_qr<SparseMatrix>
512  q (arg.sparse_matrix_value (), 0);
513  retval = ovl (q.Q (economy), q.R (economy));
514  }
515  else
516  {
517  math::sparse_qr<SparseMatrix>
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, p] = qr (a); # FIXME: not giving right dimensions.
723 %! [qe, re, pe] = qr (a, 0);
724 %! [qe2, re2, pe2] = qr (a, "econ");
725 %!
726 %! assert (q * r, a * p, sqrt (eps));
727 %! assert (qe * re, a(:, pe), sqrt (eps));
728 %! assert (qe2 * re2, a(:, pe2), sqrt (eps));
729 
730 %!test
731 %! a = [0, 2; 2, 1; 1, 2];
732 %!
733 %! [q, r] = qr (a);
734 %! [qe, re] = qr (a, 0);
735 %! [qe2, re2] = qr (a, "econ");
736 %!
737 %! assert (q * r, a, sqrt (eps));
738 %! assert (qe * re, a, sqrt (eps));
739 %! assert (qe2 * re2, a, sqrt (eps));
740 
741 %!test
742 %! a = [0, 2; 2, 1; 1, 2];
743 %!
744 %! [q, r, p] = qr (a);
745 %! [qe, re, pe] = qr (a, 0);
746 %! [qe2, re2, pe2] = qr (a, "econ");
747 %!
748 %! assert (q * r, a * p, sqrt (eps));
749 %! assert (qe * re, a(:, pe), sqrt (eps));
750 %! assert (qe2 * re2, a(:, pe2), sqrt (eps));
751 
752 %!test
753 %! a = [0, 2, 1; 2, 1, 2; 3, 1, 2];
754 %! b = [1, 3, 2; 1, 1, 0; 3, 0, 2];
755 %!
756 %! [q, r] = qr (a);
757 %! [c, re] = qr (a, b);
758 %!
759 %! assert (r, re, sqrt (eps));
760 %! assert (q'*b, c, sqrt (eps));
761 
762 %!test
763 %! a = [0, 2, i; 2, 1, 2; 3, 1, 2];
764 %! b = [1, 3, 2; 1, i, 0; 3, 0, 2];
765 %!
766 %! [q, r] = qr (a);
767 %! [c, re] = qr (a, b);
768 %!
769 %! assert (r, re, sqrt (eps));
770 %! assert (q'*b, c, sqrt (eps));
771 
772 %!test
773 %! a = [0, 2, i; 2, 1, 2; 3, 1, 2];
774 %! b = [1, 3, 2; 1, 1, 0; 3, 0, 2];
775 %!
776 %! [q, r] = qr (a);
777 %! [c, re] = qr (a, b);
778 %!
779 %! assert (r, re, sqrt (eps));
780 %! assert (q'*b, c, sqrt (eps));
781 
782 %!test
783 %! a = [0, 2, 1; 2, 1, 2; 3, 1, 2];
784 %! b = [1, 3, 2; 1, i, 0; 3, 0, 2];
785 %!
786 %! [q, r] = qr (a);
787 %! [c, re] = qr (a, b);
788 %!
789 %! assert (r, re, sqrt (eps));
790 %! assert (q'*b, c, sqrt (eps));
791 
792 ## Empty matrices
793 %!test
794 %! assert (qr (zeros (0, 0)), zeros (0, 0))
795 %! assert (qr (zeros (1, 0)), zeros (1, 0))
796 %! assert (qr (zeros (0, 1)), zeros (0, 1))
797 
798 %!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE") <*63069>
799 %! assert (qr (sparse (0, 0)), sparse (0, 0))
800 %! assert (qr (sparse (1, 0)), sparse (1, 0))
801 %! assert (qr (sparse (0, 1)), sparse (0, 1))
802 
803 %!error qr ()
804 %!error qr ([1, 2; 3, 4], 0, 2)
805 %!error <option string must be .*, not "foo"> qr (magic (3), "foo")
806 %!error <option string must be .*, not "foo"> qr (magic (3), rand (3, 1), "foo")
807 %!error <too many output arguments for dense A with B>
808 %! [q, r, p] = qr (rand (3, 2), rand (3, 1));
809 %!error <too many output arguments for dense A with B>
810 %! [q, r, p] = qr (rand (3, 2), rand (3, 1), 0);
811 
812 %!function retval = __testqr (q, r, a, p)
813 %! tol = 100* eps (class (q));
814 %! retval = 0;
815 %! if (nargin == 3)
816 %! n1 = norm (q*r - a);
817 %! n2 = norm (q'*q - eye (columns (q)));
818 %! retval = (n1 < tol && n2 < tol);
819 %! else
820 %! n1 = norm (q'*q - eye (columns (q)));
821 %! retval = (n1 < tol);
822 %! if (isvector (p))
823 %! n2 = norm (q*r - a(:,p));
824 %! retval = (retval && n2 < tol);
825 %! else
826 %! n2 = norm (q*r - a*p);
827 %! retval = (retval && n2 < tol);
828 %! endif
829 %! endif
830 %!endfunction
831 
832 %!test
833 %! t = ones (24, 1);
834 %! j = 1;
835 %!
836 %! if (false) # eliminate big matrix tests
837 %! a = rand (5000, 20);
838 %! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
839 %! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
840 %! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
841 %! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
842 %!
843 %! a = a+1i*eps;
844 %! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
845 %! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
846 %! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
847 %! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
848 %! endif
849 %!
850 %! a = [ ones(1,15); sqrt(eps)*eye(15) ];
851 %! [q,r] = qr (a); t(j++) = __testqr (q, r, a);
852 %! [q,r] = qr (a'); t(j++) = __testqr (q, r, a');
853 %! [q,r,p] = qr (a); t(j++) = __testqr (q, r, a, p);
854 %! [q,r,p] = qr (a'); t(j++) = __testqr (q, r, a', p);
855 %!
856 %! a = a+1i*eps;
857 %! [q,r] = qr (a); t(j++) = __testqr (q, r, a);
858 %! [q,r] = qr (a'); t(j++) = __testqr (q, r, a');
859 %! [q,r,p] = qr (a); t(j++) = __testqr (q, r, a, p);
860 %! [q,r,p] = qr (a'); t(j++) = __testqr (q, r, a', p);
861 %!
862 %! a = [ ones(1,15); sqrt(eps)*eye(15) ];
863 %! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
864 %! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
865 %! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
866 %! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
867 %!
868 %! a = a+1i*eps;
869 %! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
870 %! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
871 %! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
872 %! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
873 %!
874 %! a = [ 611 196 -192 407 -8 -52 -49 29
875 %! 196 899 113 -192 -71 -43 -8 -44
876 %! -192 113 899 196 61 49 8 52
877 %! 407 -192 196 611 8 44 59 -23
878 %! -8 -71 61 8 411 -599 208 208
879 %! -52 -43 49 44 -599 411 208 208
880 %! -49 -8 8 59 208 208 99 -911
881 %! 29 -44 52 -23 208 208 -911 99 ];
882 %! [q,r] = qr (a);
883 %!
884 %! assert (all (t));
885 %! if (all (t))
886 %! assert (norm (q*r - a), 0, 5000*eps);
887 %! endif
888 
889 %!test
890 %! a = single ([0, 2, 1; 2, 1, 2]);
891 %!
892 %! [q, r] = qr (a);
893 %! [qe, re] = qr (a, 0);
894 %!
895 %! assert (q * r, a, sqrt (eps ("single")));
896 %! assert (qe * re, a, sqrt (eps ("single")));
897 
898 %!test
899 %! a = single ([0, 2, 1; 2, 1, 2]);
900 %!
901 %! [q, r, p] = qr (a); # FIXME: not giving right dimensions.
902 %! [qe, re, pe] = qr (a, 0);
903 %!
904 %! assert (q * r, a * p, sqrt (eps ("single")));
905 %! assert (qe * re, a(:, pe), sqrt (eps ("single")));
906 
907 %!test
908 %! a = single ([0, 2; 2, 1; 1, 2]);
909 %!
910 %! [q, r] = qr (a);
911 %! [qe, re] = qr (a, 0);
912 %!
913 %! assert (q * r, a, sqrt (eps ("single")));
914 %! assert (qe * re, a, sqrt (eps ("single")));
915 
916 %!test
917 %! a = single ([0, 2; 2, 1; 1, 2]);
918 %!
919 %! [q, r, p] = qr (a);
920 %! [qe, re, pe] = qr (a, 0);
921 %!
922 %! assert (q * r, a * p, sqrt (eps ("single")));
923 %! assert (qe * re, a(:, pe), sqrt (eps ("single")));
924 
925 %!test
926 %! a = single ([0, 2, 1; 2, 1, 2; 3, 1, 2]);
927 %! b = single ([1, 3, 2; 1, 1, 0; 3, 0, 2]);
928 %!
929 %! [q, r] = qr (a);
930 %! [c, re] = qr (a, b);
931 %!
932 %! assert (r, re, sqrt (eps ("single")));
933 %! assert (q'*b, c, sqrt (eps ("single")));
934 
935 %!test
936 %! a = single ([0, 2, i; 2, 1, 2; 3, 1, 2]);
937 %! b = single ([1, 3, 2; 1, i, 0; 3, 0, 2]);
938 %!
939 %! [q, r] = qr (a);
940 %! [c, re] = qr (a, b);
941 %!
942 %! assert (r, re, sqrt (eps ("single")));
943 %! assert (q'*b, c, sqrt (eps ("single")));
944 
945 %!test
946 %! a = single ([0, 2, i; 2, 1, 2; 3, 1, 2]);
947 %! b = single ([1, 3, 2; 1, 1, 0; 3, 0, 2]);
948 %!
949 %! [q, r] = qr (a);
950 %! [c, re] = qr (a, b);
951 %!
952 %! assert (r, re, sqrt (eps));
953 %! assert (q'*b, c, sqrt (eps));
954 
955 %!test
956 %! a = single ([0, 2, 1; 2, 1, 2; 3, 1, 2]);
957 %! b = single ([1, 3, 2; 1, i, 0; 3, 0, 2]);
958 %!
959 %! [q, r] = qr (a);
960 %! [c, re] = qr (a, b);
961 %!
962 %! assert (r, re, sqrt (eps ("single")));
963 %! assert (q'*b, c, sqrt (eps ("single")));
964 
965 %!test
966 %! t = ones (24, 1);
967 %! j = 1;
968 %!
969 %! if (false) # eliminate big matrix tests
970 %! a = rand (5000,20);
971 %! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
972 %! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
973 %! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
974 %! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
975 %!
976 %! a = a+1i* eps ("single");
977 %! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
978 %! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
979 %! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
980 %! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
981 %! endif
982 %!
983 %! a = [ ones(1,15); sqrt(eps("single"))*eye(15) ];
984 %! [q,r] = qr (a); t(j++) = __testqr (q, r, a);
985 %! [q,r] = qr (a'); t(j++) = __testqr (q, r, a');
986 %! [q,r,p] = qr (a); t(j++) = __testqr (q, r, a, p);
987 %! [q,r,p] = qr (a'); t(j++) = __testqr (q, r, a', p);
988 %!
989 %! a = a+1i* eps ("single");
990 %! [q,r] = qr (a); t(j++) = __testqr (q, r, a);
991 %! [q,r] = qr (a'); t(j++) = __testqr (q, r, a');
992 %! [q,r,p] = qr (a); t(j++) = __testqr (q, r, a, p);
993 %! [q,r,p] = qr (a'); t(j++) = __testqr (q, r, a', p);
994 %!
995 %! a = [ ones(1,15); sqrt(eps("single"))*eye(15) ];
996 %! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
997 %! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
998 %! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
999 %! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
1000 %!
1001 %! a = a+1i* eps ("single");
1002 %! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
1003 %! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
1004 %! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
1005 %! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a',p);
1006 %!
1007 %! a = [ 611 196 -192 407 -8 -52 -49 29
1008 %! 196 899 113 -192 -71 -43 -8 -44
1009 %! -192 113 899 196 61 49 8 52
1010 %! 407 -192 196 611 8 44 59 -23
1011 %! -8 -71 61 8 411 -599 208 208
1012 %! -52 -43 49 44 -599 411 208 208
1013 %! -49 -8 8 59 208 208 99 -911
1014 %! 29 -44 52 -23 208 208 -911 99 ];
1015 %! [q,r] = qr (a);
1016 %!
1017 %! assert (all (t));
1018 %! if (all (t))
1019 %! assert (norm (q*r - a), 0, 5000 * eps ("single"));
1020 %! endif
1021 
1022 ## The deactivated tests below can't be tested till rectangular back-subs is
1023 ## implemented for sparse matrices.
1024 
1025 %!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1026 %! n = 20; d = 0.2;
1027 %! ## initialize generators to make behavior reproducible
1028 %! rand ("state", 42);
1029 %! randn ("state", 42);
1030 %! a = sprandn (n,n,d) + speye (n,n);
1031 %! r = qr (a);
1032 %! assert (r'*r, a'*a, 1e-10);
1033 
1034 %!testif HAVE_COLAMD; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1035 %! n = 20; d = 0.2;
1036 %! ## initialize generators to make behavior reproducible
1037 %! rand ("state", 42);
1038 %! randn ("state", 42);
1039 %! a = sprandn (n,n,d) + speye (n,n);
1040 %! q = symamd (a);
1041 %! a = a(q,q);
1042 %! r = qr (a);
1043 %! assert (r'*r, a'*a, 1e-10);
1044 
1045 %!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1046 %! n = 20; d = 0.2;
1047 %! ## initialize generators to make behavior reproducible
1048 %! rand ("state", 42);
1049 %! randn ("state", 42);
1050 %! a = sprandn (n,n,d) + speye (n,n);
1051 %! [c,r] = qr (a, ones (n,1));
1052 %! assert (r\c, full (a)\ones (n,1), 10e-10);
1053 
1054 %!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1055 %! n = 20; d = 0.2;
1056 %! ## initialize generators to make behavior reproducible
1057 %! rand ("state", 42);
1058 %! randn ("state", 42);
1059 %! a = sprandn (n,n,d) + speye (n,n);
1060 %! b = randn (n,2);
1061 %! [c,r] = qr (a, b);
1062 %! assert (r\c, full (a)\b, 10e-10);
1063 
1064 ## Test under-determined systems!!
1065 %!#testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1066 %! n = 20; d = 0.2;
1067 %! ## initialize generators to make behavior reproducible
1068 %! rand ("state", 42);
1069 %! randn ("state", 42);
1070 %! a = sprandn (n,n+1,d) + speye (n,n+1);
1071 %! b = randn (n,2);
1072 %! [c,r] = qr (a, b);
1073 %! assert (r\c, full (a)\b, 10e-10);
1074 
1075 %!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1076 %! n = 20; d = 0.2;
1077 %! ## initialize generators to make behavior reproducible
1078 %! rand ("state", 42);
1079 %! randn ("state", 42);
1080 %! a = 1i* sprandn (n,n,d) + speye (n,n);
1081 %! r = qr (a);
1082 %! assert (r'*r, a'*a, 1e-10);
1083 
1084 %!testif HAVE_COLAMD; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1085 %! n = 20; d = 0.2;
1086 %! ## initialize generators to make behavior reproducible
1087 %! rand ("state", 42);
1088 %! randn ("state", 42);
1089 %! a = 1i* sprandn (n,n,d) + speye (n,n);
1090 %! q = symamd (a);
1091 %! a = a(q,q);
1092 %! r = qr (a);
1093 %! assert (r'*r, a'*a, 1e-10);
1094 
1095 %!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1096 %! n = 20; d = 0.2;
1097 %! ## initialize generators to make behavior reproducible
1098 %! rand ("state", 42);
1099 %! randn ("state", 42);
1100 %! a = 1i* sprandn (n,n,d) + speye (n,n);
1101 %! [c,r] = qr (a, ones (n,1));
1102 %! assert (r\c, full (a)\ones (n,1), 10e-10);
1103 
1104 %!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1105 %! n = 20; d = 0.2;
1106 %! ## initialize generators to make behavior reproducible
1107 %! rand ("state", 42);
1108 %! randn ("state", 42);
1109 %! a = 1i* sprandn (n,n,d) + speye (n,n);
1110 %! b = randn (n,2);
1111 %! [c,r] = qr (a, b);
1112 %! assert (r\c, full (a)\b, 10e-10);
1113 
1114 ## Test under-determined systems!!
1115 %!#testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1116 %! n = 20; d = 0.2;
1117 %! ## initialize generators to make behavior reproducible
1118 %! rand ("state", 42);
1119 %! randn ("state", 42);
1120 %! a = 1i* sprandn (n,n+1,d) + speye (n,n+1);
1121 %! b = randn (n, 2);
1122 %! [c, r] = qr (a, b);
1123 %! assert (r\c, full (a)\b, 10e-10);
1124 
1125 %!testif HAVE_SPQR, HAVE_CHOLMOD
1126 %! n = 12; m = 20; d = 0.2;
1127 %! ## initialize generators to make behavior reproducible
1128 %! rand ("state", 42);
1129 %! randn ("state", 42);
1130 %! a = sprandn (m, n, d);
1131 %! b = randn (m, 2);
1132 %! [c, r] = qr (a, b);
1133 %! assert (r\c, full (a)\b, 10e-10);
1134 
1135 %!testif HAVE_SPQR, HAVE_CHOLMOD
1136 %! n = 12; m = 20; d = 0.2;
1137 %! ## initialize generators to make behavior reproducible
1138 %! rand ("state", 42);
1139 %! randn ("state", 42);
1140 %! a = sprandn (m, n, d);
1141 %! b = sprandn (m, 2, d);
1142 %! [c, r] = qr (a, b, 0);
1143 %! [c2, r2] = qr (full (a), full (b), 0);
1144 %! assert (r\c, r2\c2, 10e-10);
1145 
1146 %!testif HAVE_SPQR, HAVE_CHOLMOD
1147 %! n = 12; m = 20; d = 0.2;
1148 %! ## initialize generators to make behavior reproducible
1149 %! rand ("state", 42);
1150 %! randn ("state", 42);
1151 %! a = sprandn (m, n, d);
1152 %! b = randn (m, 2);
1153 %! [c, r, p] = qr (a, b, "matrix");
1154 %! x = p * (r\c);
1155 %! [c2, r2] = qr (full (a), b);
1156 %! x2 = r2 \ c2;
1157 %! assert (x, x2, 10e-10);
1158 
1159 %!testif HAVE_SPQR, HAVE_CHOLMOD
1160 %! n = 12; m = 20; d = 0.2;
1161 %! ## initialize generators to make behavior reproducible
1162 %! rand ("state", 42);
1163 %! randn ("state", 42);
1164 %! a = sprandn (m, n, d);
1165 %! [q, r, p] = qr (a, "matrix");
1166 %! assert (q * r, a * p, 10e-10);
1167 
1168 %!testif HAVE_SPQR, HAVE_CHOLMOD
1169 %! n = 12; m = 20; d = 0.2;
1170 %! ## initialize generators to make behavior reproducible
1171 %! rand ("state", 42);
1172 %! randn ("state", 42);
1173 %! a = sprandn (m, n, d);
1174 %! b = randn (m, 2);
1175 %! x = qr (a, b);
1176 %! [c2, r2] = qr (full (a), b);
1177 %! assert (x, r2\c2, 10e-10);
1178 
1179 %!testif HAVE_SPQR, HAVE_CHOLMOD
1180 %! n = 12; m = 20; d = 0.2;
1181 %! ## initialize generators to make behavior reproducible
1182 %! rand ("state", 42);
1183 %! randn ("state", 42);
1184 %! a = sprandn (m, n, d);
1185 %! b = i * randn (m, 2);
1186 %! x = qr (a, b);
1187 %! [c2, r2] = qr (full (a), b);
1188 %! assert (x, r2\c2, 10e-10);
1189 
1190 %!#testif HAVE_SPQR, HAVE_CHOLMOD
1191 %! n = 12; m = 20; d = 0.2;
1192 %! ## initialize generators to make behavior reproducible
1193 %! rand ("state", 42);
1194 %! randn ("state", 42);
1195 %! a = sprandn (m, n, d);
1196 %! b = i * randn (m, 2);
1197 %! [c, r] = qr (a, b);
1198 %! [c2, r2] = qr (full (a), b);
1199 %! assert (r\c, r2\c2, 10e-10);
1200 
1201 %!testif HAVE_SPQR, HAVE_CHOLMOD
1202 %! n = 12; m = 20; d = 0.2;
1203 %! ## initialize generators to make behavior reproducible
1204 %! rand ("state", 42);
1205 %! randn ("state", 42);
1206 %! a = sprandn (m, n, d);
1207 %! b = i * randn (m, 2);
1208 %! [c, r, p] = qr (a, b, "matrix");
1209 %! x = p * (r\c);
1210 %! [c2, r2] = qr (full (a), b);
1211 %! x2 = r2 \ c2;
1212 %! assert (x, x2, 10e-10);
1213 
1214 %!testif HAVE_SPQR, HAVE_CHOLMOD
1215 %! n = 12; m = 20; d = 0.2;
1216 %! ## initialize generators to make behavior reproducible
1217 %! rand ("state", 42);
1218 %! randn ("state", 42);
1219 %! a = i * sprandn (m, n, d);
1220 %! b = sprandn (m, 2, d);
1221 %! [c, r] = qr (a, b, 0);
1222 %! [c2, r2] = qr (full (a), full (b), 0);
1223 %! assert (r\c, r2\c2, 10e-10);
1224 
1225 %!testif HAVE_SPQR, HAVE_CHOLMOD
1226 %! n = 12; m = 20; d = 0.2;
1227 %! ## initialize generators to make behavior reproducible
1228 %! rand ("state", 42);
1229 %! randn ("state", 42);
1230 %! a = i * sprandn (m, n, d);
1231 %! b = randn (m, 2);
1232 %! [c, r, p] = qr (a, b, "matrix");
1233 %! x = p * (r\c);
1234 %! [c2, r2] = qr (full (a), b);
1235 %! x2 = r2 \ c2;
1236 %! assert(x, x2, 10e-10);
1237 
1238 %!testif HAVE_SPQR, HAVE_CHOLMOD
1239 %! n = 12; m = 20; d = 0.2;
1240 %! ## initialize generators to make behavior reproducible
1241 %! rand ("state", 42);
1242 %! randn ("state", 42);
1243 %! a = i * sprandn (m, n, d);
1244 %! [q, r, p] = qr (a, "matrix");
1245 %! assert(q * r, a * p, 10e-10);
1246 
1247 %!testif HAVE_SPQR, HAVE_CHOLMOD
1248 %! n = 12; m = 20; d = 0.2;
1249 %! ## initialize generators to make behavior reproducible
1250 %! rand ("state", 42);
1251 %! randn ("state", 42);
1252 %! a = i * sprandn (m, n, d);
1253 %! b = randn (m, 2);
1254 %! x = qr (a, b);
1255 %! [c2, r2] = qr (full (a), b);
1256 %! assert (x, r2\c2, 10e-10);
1257 
1258 %!testif HAVE_SPQR, HAVE_CHOLMOD
1259 %! a = sparse (5, 6);
1260 %! a(3,1) = 0.8;
1261 %! a(2,2) = 1.4;
1262 %! a(1,6) = -0.5;
1263 %! r = qr (a);
1264 %! assert (r'*r, a'*a, 10e-10);
1265 */
1266 
1267 static
1268 bool
1269 check_qr_dims (const octave_value& q, const octave_value& r,
1270  bool allow_ecf = false)
1271 {
1272  octave_idx_type m = q.rows ();
1273  octave_idx_type k = r.rows ();
1274  octave_idx_type n = r.columns ();
1275  return ((q.ndims () == 2 && r.ndims () == 2 && k == q.columns ())
1276  && (m == k || (allow_ecf && k == n && k < m)));
1277 }
1278 
1279 static
1280 bool
1281 check_index (const octave_value& i, bool vector_allowed = false)
1282 {
1283  return ((i.isreal () || i.isinteger ())
1284  && (i.is_scalar_type () || vector_allowed));
1285 }
1286 
1287 DEFUN (qrupdate, args, ,
1288  doc: /* -*- texinfo -*-
1289 @deftypefn {} {[@var{Q1}, @var{R1}] =} qrupdate (@var{Q}, @var{R}, @var{u}, @var{v})
1290 Update a QR factorization given update vectors or matrices.
1291 
1292 Given a QR@tie{}factorization of a real or complex matrix
1293 @w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and
1294 @var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization of
1295 @w{@var{A} + @var{u}*@var{v}'}, where @var{u} and @var{v} are column vectors
1296 (rank-1 update) or matrices with equal number of columns
1297 (rank-k update). Notice that the latter case is done as a sequence of
1298 rank-1 updates; thus, for k large enough, it will be both faster and more
1299 accurate to recompute the factorization from scratch.
1300 
1301 The QR@tie{}factorization supplied may be either full (Q is square) or
1302 economized (R is square).
1303 
1304 @seealso{qr, qrinsert, qrdelete, qrshift}
1305 @end deftypefn */)
1306 {
1307  octave_value_list retval;
1308 
1309  if (args.length () != 4)
1310  print_usage ();
1311 
1312  octave_value argq = args(0);
1313  octave_value argr = args(1);
1314  octave_value argu = args(2);
1315  octave_value argv = args(3);
1316 
1317  if (! argq.isnumeric () || ! argr.isnumeric ()
1318  || ! argu.isnumeric () || ! argv.isnumeric ())
1319  print_usage ();
1320 
1321  if (! check_qr_dims (argq, argr, true))
1322  error ("qrupdate: Q and R dimensions don't match");
1323 
1324  if (argq.isreal () && argr.isreal () && argu.isreal ()
1325  && argv.isreal ())
1326  {
1327  // all real case
1328  if (argq.is_single_type () || argr.is_single_type ()
1329  || argu.is_single_type () || argv.is_single_type ())
1330  {
1331  FloatMatrix Q = argq.float_matrix_value ();
1332  FloatMatrix R = argr.float_matrix_value ();
1333  FloatMatrix u = argu.float_matrix_value ();
1334  FloatMatrix v = argv.float_matrix_value ();
1335 
1336  math::qr<FloatMatrix> fact (Q, R);
1337  fact.update (u, v);
1338 
1339  retval = ovl (fact.Q (), get_qr_r (fact));
1340  }
1341  else
1342  {
1343  Matrix Q = argq.matrix_value ();
1344  Matrix R = argr.matrix_value ();
1345  Matrix u = argu.matrix_value ();
1346  Matrix v = argv.matrix_value ();
1347 
1348  math::qr<Matrix> fact (Q, R);
1349  fact.update (u, v);
1350 
1351  retval = ovl (fact.Q (), get_qr_r (fact));
1352  }
1353  }
1354  else
1355  {
1356  // complex case
1357  if (argq.is_single_type () || argr.is_single_type ()
1358  || argu.is_single_type () || argv.is_single_type ())
1359  {
1364 
1365  math::qr<FloatComplexMatrix> fact (Q, R);
1366  fact.update (u, v);
1367 
1368  retval = ovl (fact.Q (), get_qr_r (fact));
1369  }
1370  else
1371  {
1373  ComplexMatrix R = argr.complex_matrix_value ();
1374  ComplexMatrix u = argu.complex_matrix_value ();
1375  ComplexMatrix v = argv.complex_matrix_value ();
1376 
1377  math::qr<ComplexMatrix> fact (Q, R);
1378  fact.update (u, v);
1379 
1380  retval = ovl (fact.Q (), get_qr_r (fact));
1381  }
1382  }
1383 
1384  return retval;
1385 }
1386 
1387 /*
1388 %!shared A, u, v, Ac, uc, vc
1389 %! A = [0.091364 0.613038 0.999083;
1390 %! 0.594638 0.425302 0.603537;
1391 %! 0.383594 0.291238 0.085574;
1392 %! 0.265712 0.268003 0.238409;
1393 %! 0.669966 0.743851 0.445057 ];
1394 %!
1395 %! u = [0.85082;
1396 %! 0.76426;
1397 %! 0.42883;
1398 %! 0.53010;
1399 %! 0.80683 ];
1400 %!
1401 %! v = [0.98810;
1402 %! 0.24295;
1403 %! 0.43167 ];
1404 %!
1405 %! Ac = [0.620405 + 0.956953i 0.480013 + 0.048806i 0.402627 + 0.338171i;
1406 %! 0.589077 + 0.658457i 0.013205 + 0.279323i 0.229284 + 0.721929i;
1407 %! 0.092758 + 0.345687i 0.928679 + 0.241052i 0.764536 + 0.832406i;
1408 %! 0.912098 + 0.721024i 0.049018 + 0.269452i 0.730029 + 0.796517i;
1409 %! 0.112849 + 0.603871i 0.486352 + 0.142337i 0.355646 + 0.151496i ];
1410 %!
1411 %! uc = [0.20351 + 0.05401i;
1412 %! 0.13141 + 0.43708i;
1413 %! 0.29808 + 0.08789i;
1414 %! 0.69821 + 0.38844i;
1415 %! 0.74871 + 0.25821i ];
1416 %!
1417 %! vc = [0.85839 + 0.29468i;
1418 %! 0.20820 + 0.93090i;
1419 %! 0.86184 + 0.34689i ];
1420 %!
1421 
1422 %!test
1423 %! [Q,R] = qr (A);
1424 %! [Q,R] = qrupdate (Q, R, u, v);
1425 %! assert (norm (vec (Q'*Q - eye (5)), Inf), 0, 1e1*eps);
1426 %! assert (norm (vec (triu (R)-R), Inf), 0);
1427 %! assert (norm (vec (Q*R - A - u*v'), Inf), 0, norm (A)*1e1*eps);
1428 %!
1429 %!test
1430 %! [Q,R] = qr (Ac);
1431 %! [Q,R] = qrupdate (Q, R, uc, vc);
1432 %! assert (norm (vec (Q'*Q - eye (5)), Inf), 0, 1e1*eps);
1433 %! assert (norm (vec (triu (R)-R), Inf), 0);
1434 %! assert (norm (vec (Q*R - Ac - uc*vc'), Inf), 0, norm (Ac)*1e1*eps);
1435 
1436 %!test
1437 %! [Q,R] = qr (single (A));
1438 %! [Q,R] = qrupdate (Q, R, single (u), single (v));
1439 %! assert (norm (vec (Q'*Q - eye (5,"single")), Inf), single (0), ...
1440 %! 1e1 * eps ("single"));
1441 %! assert (norm (vec (triu (R)-R), Inf), single (0));
1442 %! assert (norm (vec (Q*R - single (A) - single (u)* single (v)'), Inf), ...
1443 %! single (0), norm (single (A))*1e1 * eps ("single"));
1444 %!
1445 %!test
1446 %! [Q,R] = qr (single (Ac));
1447 %! [Q,R] = qrupdate (Q, R, single (uc), single (vc));
1448 %! assert (norm (vec (Q'*Q - eye (5,"single")), Inf), single (0), ...
1449 %! 1e1 * eps ("single"));
1450 %! assert (norm (vec (triu (R)-R), Inf), single (0));
1451 %! assert (norm (vec (Q*R - single (Ac) - single (uc)* single (vc)'), Inf), ...
1452 %! single (0), norm (single (Ac))*1e1 * eps ("single"));
1453 */
1454 
1455 DEFUN (qrinsert, args, ,
1456  doc: /* -*- texinfo -*-
1457 @deftypefn {} {[@var{Q1}, @var{R1}] =} qrinsert (@var{Q}, @var{R}, @var{j}, @var{x}, @var{orient})
1458 Update a QR factorization given a row or column to insert in the original
1459 factored matrix.
1460 
1461 
1462 Given a QR@tie{}factorization of a real or complex matrix
1463 @w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and
1464 @var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization of
1465 @w{[A(:,1:j-1) x A(:,j:n)]}, where @var{u} is a column vector to be inserted
1466 into @var{A} (if @var{orient} is @qcode{"col"}), or the
1467 QR@tie{}factorization of @w{[A(1:j-1,:);x;A(:,j:n)]}, where @var{x} is a row
1468 vector to be inserted into @var{A} (if @var{orient} is @qcode{"row"}).
1469 
1470 The default value of @var{orient} is @qcode{"col"}. If @var{orient} is
1471 @qcode{"col"}, @var{u} may be a matrix and @var{j} an index vector
1472 resulting in the QR@tie{}factorization of a matrix @var{B} such that
1473 @w{B(:,@var{j})} gives @var{u} and @w{B(:,@var{j}) = []} gives @var{A}.
1474 Notice that the latter case is done as a sequence of k insertions;
1475 thus, for k large enough, it will be both faster and more accurate to
1476 recompute the factorization from scratch.
1477 
1478 If @var{orient} is @qcode{"col"}, the QR@tie{}factorization supplied may
1479 be either full (Q is square) or economized (R is square).
1480 
1481 If @var{orient} is @qcode{"row"}, full factorization is needed.
1482 @seealso{qr, qrupdate, qrdelete, qrshift}
1483 @end deftypefn */)
1484 {
1485  int nargin = args.length ();
1486 
1487  if (nargin < 4 || nargin > 5)
1488  print_usage ();
1489 
1490  octave_value argq = args(0);
1491  octave_value argr = args(1);
1492  octave_value argj = args(2);
1493  octave_value argx = args(3);
1494 
1495  if (! argq.isnumeric () || ! argr.isnumeric ()
1496  || ! argx.isnumeric ()
1497  || (nargin > 4 && ! args(4).is_string ()))
1498  print_usage ();
1499 
1500  std::string orient = (nargin < 5) ? "col" : args(4).string_value ();
1501  bool col = (orient == "col");
1502 
1503  if (! col && orient != "row")
1504  error (R"(qrinsert: ORIENT must be "col" or "row")");
1505 
1506  if (! check_qr_dims (argq, argr, col) || (! col && argx.rows () != 1))
1507  error ("qrinsert: dimension mismatch");
1508 
1509  if (! check_index (argj, col))
1510  error ("qrinsert: invalid index J");
1511 
1512  octave_value_list retval;
1513 
1515 
1516  octave_idx_type one = 1;
1517 
1518  if (argq.isreal () && argr.isreal () && argx.isreal ())
1519  {
1520  // real case
1521  if (argq.is_single_type () || argr.is_single_type ()
1522  || argx.is_single_type ())
1523  {
1524  FloatMatrix Q = argq.float_matrix_value ();
1525  FloatMatrix R = argr.float_matrix_value ();
1526  FloatMatrix x = argx.float_matrix_value ();
1527 
1528  math::qr<FloatMatrix> fact (Q, R);
1529 
1530  if (col)
1531  fact.insert_col (x, j-one);
1532  else
1533  fact.insert_row (x.row (0), j(0)-one);
1534 
1535  retval = ovl (fact.Q (), get_qr_r (fact));
1536  }
1537  else
1538  {
1539  Matrix Q = argq.matrix_value ();
1540  Matrix R = argr.matrix_value ();
1541  Matrix x = argx.matrix_value ();
1542 
1543  math::qr<Matrix> fact (Q, R);
1544 
1545  if (col)
1546  fact.insert_col (x, j-one);
1547  else
1548  fact.insert_row (x.row (0), j(0)-one);
1549 
1550  retval = ovl (fact.Q (), get_qr_r (fact));
1551  }
1552  }
1553  else
1554  {
1555  // complex case
1556  if (argq.is_single_type () || argr.is_single_type ()
1557  || argx.is_single_type ())
1558  {
1562 
1563  math::qr<FloatComplexMatrix> fact (Q, R);
1564 
1565  if (col)
1566  fact.insert_col (x, j-one);
1567  else
1568  fact.insert_row (x.row (0), j(0)-one);
1569 
1570  retval = ovl (fact.Q (), get_qr_r (fact));
1571  }
1572  else
1573  {
1575  ComplexMatrix R = argr.complex_matrix_value ();
1577 
1578  math::qr<ComplexMatrix> fact (Q, R);
1579 
1580  if (col)
1581  fact.insert_col (x, j-one);
1582  else
1583  fact.insert_row (x.row (0), j(0)-one);
1584 
1585  retval = ovl (fact.Q (), get_qr_r (fact));
1586  }
1587  }
1588 
1589  return retval;
1590 }
1591 
1592 /*
1593 %!test
1594 %! [Q,R] = qr (A);
1595 %! [Q,R] = qrinsert (Q, R, 3, u);
1596 %! assert (norm (vec (Q'*Q - eye (5)), Inf), 0, 1e1*eps);
1597 %! assert (norm (vec (triu (R) - R), Inf), 0);
1598 %! assert (norm (vec (Q*R - [A(:,1:2) u A(:,3)]), Inf), 0, norm (A)*1e1*eps);
1599 %!test
1600 %! [Q,R] = qr (Ac);
1601 %! [Q,R] = qrinsert (Q, R, 3, uc);
1602 %! assert (norm (vec (Q'*Q - eye (5)), Inf), 0, 1e1*eps);
1603 %! assert (norm (vec (triu (R) - R), Inf), 0);
1604 %! assert (norm (vec (Q*R - [Ac(:,1:2) uc Ac(:,3)]), Inf), 0, ...
1605 %! 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), 0, 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), 0, 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), 0, 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), 0, ...
1622 %! norm (Ac) * 1e1 * eps);
1623 
1624 %!test
1625 %! [Q,R] = qr (single (A));
1626 %! [Q,R] = qrinsert (Q, R, 3, single (u));
1627 %! assert (norm (vec (Q'*Q - eye (5,"single")), Inf), single (0), ...
1628 %! 1e1 * eps ("single"));
1629 %! assert (norm (vec (triu (R) - R), Inf), single (0));
1630 %! assert (norm (vec (Q*R - single ([A(:,1:2) u A(:,3)])), Inf), ...
1631 %! single (0), norm (single (A))*1e1 * eps ("single"));
1632 %!test
1633 %! [Q,R] = qr (single (Ac));
1634 %! [Q,R] = qrinsert (Q, R, 3, single (uc));
1635 %! assert (norm (vec (Q'*Q - eye (5,"single")), Inf), single (0), ...
1636 %! 1e1 * eps ("single"));
1637 %! assert (norm (vec (triu (R) - R), Inf), single (0));
1638 %! assert (norm (vec (Q*R - single ([Ac(:,1:2) uc Ac(:,3)])), Inf), ...
1639 %! single (0), norm (single (Ac))*1e1 * eps ("single"));
1640 %!test
1641 %! x = single ([0.85082 0.76426 0.42883 ]);
1642 %!
1643 %! [Q,R] = qr (single (A));
1644 %! [Q,R] = qrinsert (Q, R, 3, x, "row");
1645 %! assert (norm (vec (Q'*Q - eye (6,"single")), Inf), single (0), ...
1646 %! 1e1 * eps ("single"));
1647 %! assert (norm (vec (triu (R) - R), Inf), single (0));
1648 %! assert (norm (vec (Q*R - single ([A(1:2,:);x;A(3:5,:)])), Inf), ...
1649 %! single (0), norm (single (A))*1e1 * eps ("single"));
1650 %!test
1651 %! x = single ([0.20351 + 0.05401i 0.13141 + 0.43708i 0.29808 + 0.08789i ]);
1652 %!
1653 %! [Q,R] = qr (single (Ac));
1654 %! [Q,R] = qrinsert (Q, R, 3, x, "row");
1655 %! assert (norm (vec (Q'*Q - eye (6,"single")), Inf), single (0), ...
1656 %! 1e1 * eps ("single"));
1657 %! assert (norm (vec (triu (R) - R), Inf), single (0));
1658 %! assert (norm (vec (Q*R - single ([Ac(1:2,:);x;Ac(3:5,:)])), Inf), ...
1659 %! single (0), norm (single (Ac))*1e1 * eps ("single"));
1660 */
1661 
1662 DEFUN (qrdelete, args, ,
1663  doc: /* -*- texinfo -*-
1664 @deftypefn {} {[@var{Q1}, @var{R1}] =} qrdelete (@var{Q}, @var{R}, @var{j}, @var{orient})
1665 Update a QR factorization given a row or column to delete from the original
1666 factored matrix.
1667 
1668 Given a QR@tie{}factorization of a real or complex matrix
1669 @w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and
1670 @var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization of
1671 @w{[A(:,1:j-1), U, A(:,j:n)]},
1672 where @var{u} is a column vector to be inserted into @var{A}
1673 (if @var{orient} is @qcode{"col"}),
1674 or the QR@tie{}factorization of @w{[A(1:j-1,:);X;A(:,j:n)]},
1675 where @var{x} is a row @var{orient} is @qcode{"row"}).
1676 The default value of @var{orient} is @qcode{"col"}.
1677 
1678 If @var{orient} is @qcode{"col"}, @var{j} may be an index vector
1679 resulting in the QR@tie{}factorization of a matrix @var{B} such that
1680 @w{A(:,@var{j}) = []} gives @var{B}. Notice that the latter case is done as
1681 a sequence of k deletions; thus, for k large enough, it will be both faster
1682 and more accurate to recompute the factorization from scratch.
1683 
1684 If @var{orient} is @qcode{"col"}, the QR@tie{}factorization supplied may
1685 be either full (Q is square) or economized (R is square).
1686 
1687 If @var{orient} is @qcode{"row"}, full factorization is needed.
1688 @seealso{qr, qrupdate, qrinsert, qrshift}
1689 @end deftypefn */)
1690 {
1691  int nargin = args.length ();
1692 
1693  if (nargin < 3 || nargin > 4)
1694  print_usage ();
1695 
1696  octave_value argq = args(0);
1697  octave_value argr = args(1);
1698  octave_value argj = args(2);
1699 
1700  if (! argq.isnumeric () || ! argr.isnumeric ()
1701  || (nargin > 3 && ! args(3).is_string ()))
1702  print_usage ();
1703 
1704  std::string orient = (nargin < 4) ? "col" : args(3).string_value ();
1705  bool col = orient == "col";
1706 
1707  if (! col && orient != "row")
1708  error (R"(qrdelete: ORIENT must be "col" or "row")");
1709 
1710  if (! check_qr_dims (argq, argr, col))
1711  error ("qrdelete: dimension mismatch");
1712 
1714  if (! check_index (argj, col))
1715  error ("qrdelete: invalid index J");
1716 
1717  octave_value_list retval;
1718 
1719  octave_idx_type one = 1;
1720 
1721  if (argq.isreal () && argr.isreal ())
1722  {
1723  // real case
1724  if (argq.is_single_type () || argr.is_single_type ())
1725  {
1726  FloatMatrix Q = argq.float_matrix_value ();
1727  FloatMatrix R = argr.float_matrix_value ();
1728 
1729  math::qr<FloatMatrix> fact (Q, R);
1730 
1731  if (col)
1732  fact.delete_col (j-one);
1733  else
1734  fact.delete_row (j(0)-one);
1735 
1736  retval = ovl (fact.Q (), get_qr_r (fact));
1737  }
1738  else
1739  {
1740  Matrix Q = argq.matrix_value ();
1741  Matrix R = argr.matrix_value ();
1742 
1743  math::qr<Matrix> fact (Q, R);
1744 
1745  if (col)
1746  fact.delete_col (j-one);
1747  else
1748  fact.delete_row (j(0)-one);
1749 
1750  retval = ovl (fact.Q (), get_qr_r (fact));
1751  }
1752  }
1753  else
1754  {
1755  // complex case
1756  if (argq.is_single_type () || argr.is_single_type ())
1757  {
1760 
1761  math::qr<FloatComplexMatrix> fact (Q, R);
1762 
1763  if (col)
1764  fact.delete_col (j-one);
1765  else
1766  fact.delete_row (j(0)-one);
1767 
1768  retval = ovl (fact.Q (), get_qr_r (fact));
1769  }
1770  else
1771  {
1773  ComplexMatrix R = argr.complex_matrix_value ();
1774 
1775  math::qr<ComplexMatrix> fact (Q, R);
1776 
1777  if (col)
1778  fact.delete_col (j-one);
1779  else
1780  fact.delete_row (j(0)-one);
1781 
1782  retval = ovl (fact.Q (), get_qr_r (fact));
1783  }
1784  }
1785 
1786  return retval;
1787 }
1788 
1789 /*
1790 %!test
1791 %! AA = [0.091364 0.613038 0.027504 0.999083;
1792 %! 0.594638 0.425302 0.562834 0.603537;
1793 %! 0.383594 0.291238 0.742073 0.085574;
1794 %! 0.265712 0.268003 0.783553 0.238409;
1795 %! 0.669966 0.743851 0.457255 0.445057 ];
1796 %!
1797 %! [Q,R] = qr (AA);
1798 %! [Q,R] = qrdelete (Q, R, 3);
1799 %! assert (norm (vec (Q'*Q - eye (5)), Inf), 0, 16*eps);
1800 %! assert (norm (vec (triu (R) - R), Inf), 0);
1801 %! assert (norm (vec (Q*R - [AA(:,1:2) AA(:,4)]), Inf), 0, norm (AA)*1e1*eps);
1802 %!
1803 %!test
1804 %! AA = [0.364554 + 0.993117i 0.669818 + 0.510234i 0.426568 + 0.041337i 0.847051 + 0.233291i;
1805 %! 0.049600 + 0.242783i 0.448946 + 0.484022i 0.141155 + 0.074420i 0.446746 + 0.392706i;
1806 %! 0.581922 + 0.657416i 0.581460 + 0.030016i 0.219909 + 0.447288i 0.201144 + 0.069132i;
1807 %! 0.694986 + 0.000571i 0.682327 + 0.841712i 0.807537 + 0.166086i 0.192767 + 0.358098i;
1808 %! 0.945002 + 0.066788i 0.350492 + 0.642638i 0.579629 + 0.048102i 0.600170 + 0.636938i ] * I;
1809 %!
1810 %! [Q,R] = qr (AA);
1811 %! [Q,R] = qrdelete (Q, R, 3);
1812 %! assert (norm (vec (Q'*Q - eye (5)), Inf) < 16*eps);
1813 %! assert (norm (vec (triu (R) - R), Inf) == 0);
1814 %! assert (norm (vec (Q*R - [AA(:,1:2) AA(:,4)]), Inf) < norm (AA)*1e1*eps);
1815 %!
1816 %!test
1817 %! AA = [0.091364 0.613038 0.027504 0.999083;
1818 %! 0.594638 0.425302 0.562834 0.603537;
1819 %! 0.383594 0.291238 0.742073 0.085574;
1820 %! 0.265712 0.268003 0.783553 0.238409;
1821 %! 0.669966 0.743851 0.457255 0.445057 ];
1822 %!
1823 %! [Q,R] = qr (AA);
1824 %! [Q,R] = qrdelete (Q, R, 3, "row");
1825 %! assert (norm (vec (Q'*Q - eye (4)), Inf), 0, 1e1*eps);
1826 %! assert (norm (vec (triu (R) - R), Inf), 0);
1827 %! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf), 0, norm (AA)*1e1*eps);
1828 %!
1829 %!test
1830 %! AA = [0.364554 + 0.993117i 0.669818 + 0.510234i 0.426568 + 0.041337i 0.847051 + 0.233291i;
1831 %! 0.049600 + 0.242783i 0.448946 + 0.484022i 0.141155 + 0.074420i 0.446746 + 0.392706i;
1832 %! 0.581922 + 0.657416i 0.581460 + 0.030016i 0.219909 + 0.447288i 0.201144 + 0.069132i;
1833 %! 0.694986 + 0.000571i 0.682327 + 0.841712i 0.807537 + 0.166086i 0.192767 + 0.358098i;
1834 %! 0.945002 + 0.066788i 0.350492 + 0.642638i 0.579629 + 0.048102i 0.600170 + 0.636938i ] * I;
1835 %!
1836 %! [Q,R] = qr (AA);
1837 %! [Q,R] = qrdelete (Q, R, 3, "row");
1838 %! assert (norm (vec (Q'*Q - eye (4)), Inf), 0, 1e1*eps);
1839 %! assert (norm (vec (triu (R) - R), Inf), 0);
1840 %! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf), 0, norm (AA)*1e1*eps);
1841 
1842 %!test
1843 %! AA = single ([0.091364 0.613038 0.027504 0.999083;
1844 %! 0.594638 0.425302 0.562834 0.603537;
1845 %! 0.383594 0.291238 0.742073 0.085574;
1846 %! 0.265712 0.268003 0.783553 0.238409;
1847 %! 0.669966 0.743851 0.457255 0.445057 ]);
1848 %!
1849 %! [Q,R] = qr (AA);
1850 %! [Q,R] = qrdelete (Q, R, 3);
1851 %! assert (norm (vec (Q'*Q - eye (5,"single")), Inf), single (0), ...
1852 %! 1e1 * eps ("single"));
1853 %! assert (norm (vec (triu (R) - R), Inf), single (0));
1854 %! assert (norm (vec (Q*R - [AA(:,1:2) AA(:,4)]), Inf), single (0), ...
1855 %! norm (AA)*1e1 * eps ("single"));
1856 %!
1857 %!test
1858 %! AA = single ([0.364554 + 0.993117i 0.669818 + 0.510234i 0.426568 + 0.041337i 0.847051 + 0.233291i;
1859 %! 0.049600 + 0.242783i 0.448946 + 0.484022i 0.141155 + 0.074420i 0.446746 + 0.392706i;
1860 %! 0.581922 + 0.657416i 0.581460 + 0.030016i 0.219909 + 0.447288i 0.201144 + 0.069132i;
1861 %! 0.694986 + 0.000571i 0.682327 + 0.841712i 0.807537 + 0.166086i 0.192767 + 0.358098i;
1862 %! 0.945002 + 0.066788i 0.350492 + 0.642638i 0.579629 + 0.048102i 0.600170 + 0.636938i ]) * I;
1863 %!
1864 %! [Q,R] = qr (AA);
1865 %! [Q,R] = qrdelete (Q, R, 3);
1866 %! assert (norm (vec (Q'*Q - eye (5,"single")), Inf), single (0), ...
1867 %! 1e1 * eps ("single"));
1868 %! assert (norm (vec (triu (R) - R), Inf), single (0));
1869 %! assert (norm (vec (Q*R - [AA(:,1:2) AA(:,4)]), Inf), single (0), ...
1870 %! norm (AA)*1e1 * eps ("single"));
1871 
1872 %!test
1873 %! AA = single ([0.091364 0.613038 0.027504 0.999083;
1874 %! 0.594638 0.425302 0.562834 0.603537;
1875 %! 0.383594 0.291238 0.742073 0.085574;
1876 %! 0.265712 0.268003 0.783553 0.238409;
1877 %! 0.669966 0.743851 0.457255 0.445057 ]);
1878 %!
1879 %! [Q,R] = qr (AA);
1880 %! [Q,R] = qrdelete (Q, R, 3, "row");
1881 %! assert (norm (vec (Q'*Q - eye (4,"single")), Inf), single (0), ...
1882 %! 1.5e1 * eps ("single"));
1883 %! assert (norm (vec (triu (R) - R), Inf), single (0));
1884 %! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf), single (0), ...
1885 %! norm (AA)*1e1 * eps ("single"));
1886 %!testif HAVE_QRUPDATE
1887 %! ## Same test as above but with more precicision
1888 %! AA = single ([0.091364 0.613038 0.027504 0.999083;
1889 %! 0.594638 0.425302 0.562834 0.603537;
1890 %! 0.383594 0.291238 0.742073 0.085574;
1891 %! 0.265712 0.268003 0.783553 0.238409;
1892 %! 0.669966 0.743851 0.457255 0.445057 ]);
1893 %!
1894 %! [Q,R] = qr (AA);
1895 %! [Q,R] = qrdelete (Q, R, 3, "row");
1896 %! assert (norm (vec (Q'*Q - eye (4,"single")), Inf), single (0), ...
1897 %! 1e1* eps ("single"));
1898 %! assert (norm (vec (triu (R) - R), Inf), single (0));
1899 %! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf), single (0), ...
1900 %! norm (AA)*1e1 * eps ("single"));
1901 %!
1902 %!test
1903 %! AA = single ([0.364554 + 0.993117i 0.669818 + 0.510234i 0.426568 + 0.041337i 0.847051 + 0.233291i;
1904 %! 0.049600 + 0.242783i 0.448946 + 0.484022i 0.141155 + 0.074420i 0.446746 + 0.392706i;
1905 %! 0.581922 + 0.657416i 0.581460 + 0.030016i 0.219909 + 0.447288i 0.201144 + 0.069132i;
1906 %! 0.694986 + 0.000571i 0.682327 + 0.841712i 0.807537 + 0.166086i 0.192767 + 0.358098i;
1907 %! 0.945002 + 0.066788i 0.350492 + 0.642638i 0.579629 + 0.048102i 0.600170 + 0.636938i ]) * I;
1908 %!
1909 %! [Q,R] = qr (AA);
1910 %! [Q,R] = qrdelete (Q, R, 3, "row");
1911 %! assert (norm (vec (Q'*Q - eye (4,"single")), Inf), single (0), ...
1912 %! 1e1 * eps ("single"));
1913 %! assert (norm (vec (triu (R) - R), Inf), single (0));
1914 %! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf), single (0), ...
1915 %! norm (AA)*1e1 * eps ("single"));
1916 */
1917 
1918 DEFUN (qrshift, args, ,
1919  doc: /* -*- texinfo -*-
1920 @deftypefn {} {[@var{Q1}, @var{R1}] =} qrshift (@var{Q}, @var{R}, @var{i}, @var{j})
1921 Update a QR factorization given a range of columns to shift in the original
1922 factored matrix.
1923 
1924 Given a QR@tie{}factorization of a real or complex matrix
1925 @w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and
1926 @var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization
1927 of @w{@var{A}(:,p)}, where @w{p} is the permutation @*
1928 @code{p = [1:i-1, shift(i:j, 1), j+1:n]} if @w{@var{i} < @var{j}} @*
1929  or @*
1930 @code{p = [1:j-1, shift(j:i,-1), i+1:n]} if @w{@var{j} < @var{i}}. @*
1931 
1932 @seealso{qr, qrupdate, qrinsert, qrdelete}
1933 @end deftypefn */)
1934 {
1935  if (args.length () != 4)
1936  print_usage ();
1937 
1938  octave_value argq = args(0);
1939  octave_value argr = args(1);
1940  octave_value argi = args(2);
1941  octave_value argj = args(3);
1942 
1943  if (! argq.isnumeric () || ! argr.isnumeric ())
1944  print_usage ();
1945 
1946  if (! check_qr_dims (argq, argr, true))
1947  error ("qrshift: dimensions mismatch");
1948 
1949  octave_idx_type i = argi.idx_type_value ();
1950  octave_idx_type j = argj.idx_type_value ();
1951 
1952  if (! check_index (argi) || ! check_index (argj))
1953  error ("qrshift: invalid index I or J");
1954 
1955  octave_value_list retval;
1956 
1957  if (argq.isreal () && argr.isreal ())
1958  {
1959  // all real case
1960  if (argq.is_single_type ()
1961  && argr.is_single_type ())
1962  {
1963  FloatMatrix Q = argq.float_matrix_value ();
1964  FloatMatrix R = argr.float_matrix_value ();
1965 
1966  math::qr<FloatMatrix> fact (Q, R);
1967  fact.shift_cols (i-1, j-1);
1968 
1969  retval = ovl (fact.Q (), get_qr_r (fact));
1970  }
1971  else
1972  {
1973  Matrix Q = argq.matrix_value ();
1974  Matrix R = argr.matrix_value ();
1975 
1976  math::qr<Matrix> fact (Q, R);
1977  fact.shift_cols (i-1, j-1);
1978 
1979  retval = ovl (fact.Q (), get_qr_r (fact));
1980  }
1981  }
1982  else
1983  {
1984  // complex case
1985  if (argq.is_single_type ()
1986  && argr.is_single_type ())
1987  {
1990 
1991  math::qr<FloatComplexMatrix> fact (Q, R);
1992  fact.shift_cols (i-1, j-1);
1993 
1994  retval = ovl (fact.Q (), get_qr_r (fact));
1995  }
1996  else
1997  {
1999  ComplexMatrix R = argr.complex_matrix_value ();
2000 
2001  math::qr<ComplexMatrix> fact (Q, R);
2002  fact.shift_cols (i-1, j-1);
2003 
2004  retval = ovl (fact.Q (), get_qr_r (fact));
2005  }
2006  }
2007 
2008  return retval;
2009 }
2010 
2011 /*
2012 %!test
2013 %! AA = A.';
2014 %! i = 2; j = 4; p = [1:i-1, shift(i:j,-1), j+1:5];
2015 %!
2016 %! [Q,R] = qr (AA);
2017 %! [Q,R] = qrshift (Q, R, i, j);
2018 %! assert (norm (vec (Q'*Q - eye (3)), Inf), 0, 1e1*eps);
2019 %! assert (norm (vec (triu (R) - R), Inf), 0);
2020 %! assert (norm (vec (Q*R - AA(:,p)), Inf), 0, norm (AA)*1e1*eps);
2021 %!
2022 %! j = 2; i = 4; p = [1:j-1, shift(j:i,+1), i+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), 0, 1e1*eps);
2027 %! assert (norm (vec (triu (R) - R), Inf), 0);
2028 %! assert (norm (vec (Q*R - AA(:,p)), Inf), 0, norm (AA)*1e1*eps);
2029 %!
2030 %!test
2031 %! AA = Ac.';
2032 %! i = 2; j = 4; p = [1:i-1, shift(i:j,-1), j+1:5];
2033 %!
2034 %! [Q,R] = qr (AA);
2035 %! [Q,R] = qrshift (Q, R, i, j);
2036 %! assert (norm (vec (Q'*Q - eye (3)), Inf), 0, 1e1*eps);
2037 %! assert (norm (vec (triu (R) - R), Inf), 0);
2038 %! assert (norm (vec (Q*R - AA(:,p)), Inf), 0, norm (AA)*1e1*eps);
2039 %!
2040 %! j = 2; i = 4; p = [1:j-1, shift(j:i,+1), i+1:5];
2041 %!
2042 %! [Q,R] = qr (AA);
2043 %! [Q,R] = qrshift (Q, R, i, j);
2044 %! assert (norm (vec (Q'*Q - eye (3)), Inf), 0, 1e1*eps);
2045 %! assert (norm (vec (triu (R) - R), Inf), 0);
2046 %! assert (norm (vec (Q*R - AA(:,p)), Inf), 0, norm (AA)*1e1*eps);
2047 
2048 %!test
2049 %! AA = single (A).';
2050 %! i = 2; j = 4; p = [1:i-1, shift(i:j,-1), j+1:5];
2051 %!
2052 %! [Q,R] = qr (AA);
2053 %! [Q,R] = qrshift (Q, R, i, j);
2054 %! assert (norm (vec (Q'*Q - eye (3,"single")), Inf), single (0), ...
2055 %! 1e1 * eps ("single"));
2056 %! assert (norm (vec (triu (R) - R), Inf), single (0));
2057 %! assert (norm (vec (Q*R - AA(:,p)), Inf), single (0), ...
2058 %! norm (AA)*1e1 * eps ("single"));
2059 %!
2060 %! j = 2; i = 4; p = [1:j-1, shift(j:i,+1), i+1:5];
2061 %!
2062 %! [Q,R] = qr (AA);
2063 %! [Q,R] = qrshift (Q, R, i, j);
2064 %! assert (norm (vec (Q'*Q - eye (3,"single")), Inf), single (0), ...
2065 %! 1e1 * eps ("single"));
2066 %! assert (norm (vec (triu (R) - R), Inf), single (0));
2067 %! assert (norm (vec (Q*R - AA(:,p)), Inf), single (0), ...
2068 %! norm (AA)*1e1 * eps ("single"));
2069 %!
2070 %!test
2071 %! AA = single (Ac).';
2072 %! i = 2; j = 4; p = [1:i-1, shift(i:j,-1), j+1:5];
2073 %!
2074 %! [Q,R] = qr (AA);
2075 %! [Q,R] = qrshift (Q, R, i, j);
2076 %! assert (norm (vec (Q'*Q - eye (3,"single")), Inf), single (0), ...
2077 %! 1e1 * eps ("single"));
2078 %! assert (norm (vec (triu (R) - R), Inf), single (0));
2079 %! assert (norm (vec (Q*R - AA(:,p)), Inf), single (0), ...
2080 %! norm (AA)*1e1 * eps ("single"));
2081 %!
2082 %! j = 2; i = 4; p = [1:j-1, shift(j:i,+1), i+1:5];
2083 %!
2084 %! [Q,R] = qr (AA);
2085 %! [Q,R] = qrshift (Q, R, i, j);
2086 %! assert (norm (vec (Q'*Q - eye (3,"single")), Inf), single (0), ...
2087 %! 1e1 * eps ("single"));
2088 %! assert (norm (vec (triu (R) - R), Inf), single (0));
2089 %! assert (norm (vec (Q*R - AA(:,p)), Inf), single (0), ...
2090 %! norm (AA)*1e1 * eps ("single"));
2091 */
2092 
2093 OCTAVE_END_NAMESPACE(octave)
ComplexColumnVector conj(const ComplexColumnVector &a)
Definition: CColVector.cc:217
Definition: dMatrix.h:42
bool isinteger() const
Definition: ov.h:730
SparseMatrix sparse_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:900
int ndims() const
Definition: ov.h:551
octave_idx_type rows() const
Definition: ov.h:545
bool is_scalar_type() const
Definition: ov.h:744
bool isreal() const
Definition: ov.h:738
bool isnumeric() const
Definition: ov.h:750
ComplexMatrix complex_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:871
octave_idx_type idx_type_value(bool req_int=false, bool frc_str_conv=false) const
bool is_single_type() const
Definition: ov.h:698
Array< octave_idx_type > octave_idx_type_vector_value(bool req_int=false, bool frc_str_conv=false, bool frc_vec_conv=false) const
bool issparse() const
Definition: ov.h:753
FloatMatrix float_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:856
bool iscomplex() const
Definition: ov.h:741
octave_idx_type columns() const
Definition: ov.h:547
Matrix matrix_value(bool frc_str_conv=false) const
Definition: ov.h:853
FloatComplexMatrix float_complex_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:875
SparseComplexMatrix sparse_complex_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:904
Definition: qr.h:40
OCTAVE_BEGIN_NAMESPACE(octave) static octave_value daspk_fcn
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:988
void err_wrong_type_arg(const char *name, const char *s)
Definition: errwarn.cc:166
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:781
octave_idx_type n
Definition: mx-inlines.cc:761
T * r
Definition: mx-inlines.cc:781
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:219