GNU Octave 10.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
 
Loading...
Searching...
No Matches
qr.cc
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////
2//
3// Copyright (C) 1996-2025 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
61template <typename MT>
62static octave_value
63get_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
72template <typename T>
73static typename math::qr<T>::type
74qr_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//
110DEFUN (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{}, "econ")
119@deftypefnx {} {[@dots{}] =} qr (@dots{}, "vector")
120@deftypefnx {} {[@dots{}] =} qr (@dots{}, "matrix")
121@deftypefnx {} {[@dots{}] =} qr (@dots{}, 0)
122@cindex QR factorization
123Compute the QR@tie{}factorization of @var{A}, using standard @sc{lapack}
124subroutines.
125
126The 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
137where @var{Q} is an orthogonal matrix and @var{R} is upper triangular.
138@end ifnottex
139
140For 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
147returns
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
164which 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
175If 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
178the first return value when multiple values are requested.)
179
180If a third output @var{P} is requested, then @code{qr} calculates the permuted
181QR@tie{}factorization
182@tex
183$QR = AP$ where $Q$ is an orthogonal matrix, $R$ is upper triangular, and $P$
184is 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
193where @var{Q} is an orthogonal matrix, @var{R} is upper triangular, and
194@var{P} is a permutation matrix.
195@end ifnottex
196
197If @var{A} is dense, the permuted QR@tie{}factorization has the additional
198property that the diagonal entries of @var{R} are ordered by decreasing
199magnitude. In other words, @code{abs (diag (@var{R}))} will be ordered
200from largest to smallest.
201
202If @var{A} is sparse, @var{P} is a fill-reducing ordering of the columns
203of @var{A}. In that case, the diagonal entries of @var{R} are not ordered by
204decreasing magnitude.
205
206For 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
213returns
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
234If the input matrix @var{A} is sparse, the sparse QR@tie{}factorization
235is computed by using @sc{SPQR} or @sc{cxsparse} (e.g., if @sc{SPQR} is not
236available). Because the matrix @var{Q} is, in general, a full matrix, it is
237recommended to request only one return value @var{R}. In that case, the
238computation avoids the construction of @var{Q} and returns a sparse @var{R}
239such that @code{@var{R} = chol (@var{A}' * @var{A})}.
240
241If @var{A} is dense, an additional input matrix @var{B} is supplied, and two
242return values are requested, then @code{qr} returns @var{C}, where
243@code{@var{C} = @var{Q}' * @var{B}}. This allows the least squares
244approximation 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
253If @var{A} is a sparse @nospell{MxN} matrix and an additional matrix @var{B} is
254supplied, one or two return values are possible. If one return value @var{X}
255is 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
257approximation @w{of @code{@var{A} \ @var{B}}}. If two return values are
258requested, @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
260parameter should be preferred because it uses less memory and can handle
261rank-deficient matrices better.
262
263If the final argument is the string @qcode{"vector"} then @var{P} is a
264permutation vector (of the columns of @var{A}) instead of a permutation
265matrix. In this case, the defining relationship is:
266
267@example
268@var{Q} * @var{R} = @var{A}(:, @var{P})
269@end example
270
271The default, however, is to return a permutation matrix and this may be
272explicitly specified by using a final argument of @qcode{"matrix"}.
273
274When the optional argument is the string @qcode{"econ"}, an economy
275factorization is returned. If the original matrix @var{A} has size
276@nospell{MxN} and M > N, then the economy factorization will calculate just N
277rows 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
279factorizations.
280
281If the optional argument is the numeric value 0 then @code{qr} acts as if
282the @qcode{"econ"} and @qcode{"vector"} arguments were both given.
283@strong{Warning:} This syntax is accepted, but no longer recommended and may
284be removed in the future. Use @qcode{"econ"} instead.
285
286Background: The QR factorization has applications in the solution of least
287squares problems
288@tex
289$$
290\min_x \left\Vert A x - b \right\Vert_2
291$$
292@end tex
293@ifnottex
294
295@example
296min norm (A*x - b)
297@end example
298
299@end ifnottex
300for overdetermined systems of equations (i.e.,
301@tex
302$A$
303@end tex
304@ifnottex
305@var{A}
306@end ifnottex
307is a tall, thin matrix).
308
309The permuted QR@tie{}factorization
310@code{[@var{Q}, @var{R}, @var{P}] = qr (@var{A})} allows the construction of an
311orthogonal basis of @code{span (A)}.
312
313@seealso{chol, hess, lu, qz, schur, svd, qrupdate, qrinsert, qrdelete, qrshift}
314@end deftypefn */)
315{
316 int nargin = args.length ();
317
318 if (nargin < 1 || nargin > 3)
319 print_usage ();
320
321 if (nargout > 3)
322 error ("qr: too many output arguments");
323
324 octave_value_list retval;
325
326 octave_value arg = args(0);
327
328 bool economy = false;
329 bool is_cmplx = false;
330 bool have_b = false;
331 bool vector_p = false;
332
333 if (arg.iscomplex ())
334 is_cmplx = true;
335 if (nargin > 1)
336 {
337 have_b = true;
338 if (args(nargin-1).is_scalar_type ())
339 {
340 int val = args(nargin-1).int_value ();
341 if (val == 0)
342 {
343 economy = true;
344 vector_p = true;
345 have_b = (nargin > 2);
346 }
347 else if (nargin == 3) // argument 3 should be 0 or a string
348 print_usage ();
349 }
350 else if (args(nargin-1).is_string ())
351 {
352 std::string str = args(nargin-1).string_value ();
353 if (str == "econ")
354 economy = true;
355 else if (str == "vector")
356 vector_p = true;
357 else if (str != "matrix")
358 error ("qr: option string must be 'econ' or 'matrix' or 'vector', not \"%s\"", str.c_str ());
359 have_b = (nargin > 2);
360 }
361 else if (! args(nargin-1).isnumeric ())
362 err_wrong_type_arg ("qr", args(nargin-1));
363 else if (nargin == 3) // should be caught by is_scalar_type or is_string
364 print_usage ();
365
366 if (have_b && ! args(1).isnumeric ())
367 print_usage ();
368
369 if (have_b && args(1).iscomplex ())
370 is_cmplx = true;
371 }
372
373 if (arg.issparse ())
374 {
375 if (is_cmplx)
376 {
377 if (have_b && nargout == 1)
378 {
379 octave_idx_type info;
380
381 if (! args(1).issparse () && args(1).iscomplex ())
382 retval = ovl
383 (math::sparse_qr<SparseComplexMatrix>::solve
386 args(1).complex_matrix_value (), info));
387 else if (args(1).issparse () && args(1).iscomplex ())
388 retval = ovl
389 (math::sparse_qr<SparseComplexMatrix>::solve
392 args(1).sparse_complex_matrix_value (), info));
393 else if (! args(1).issparse () && ! args(1).iscomplex ())
394 retval = ovl
395 (math::sparse_qr<SparseComplexMatrix>::solve
398 args(1).matrix_value (), info));
399 else if (args(1).issparse () && ! args(1).iscomplex ())
400 retval = ovl
401 (math::sparse_qr<SparseComplexMatrix>::solve
404 args(1).sparse_matrix_value (), info));
405 }
406 else if (have_b && nargout == 2)
407 {
408 math::sparse_qr<SparseComplexMatrix>
409 q (arg.sparse_complex_matrix_value (), 0);
410 retval = ovl (q.C (args(1).complex_matrix_value (), economy),
411 q.R (economy));
412 }
413 else if (have_b && nargout == 3)
414 {
415 math::sparse_qr<SparseComplexMatrix>
417 if (vector_p)
418 retval = ovl (q.C (args(1).complex_matrix_value (), economy),
419 q.R (economy), q.E ().transpose ());
420 else
421 retval = ovl (q.C (args(1).complex_matrix_value (), economy),
422 q.R (economy), q.E_MAT ());
423 }
424 else
425 {
426 if (nargout > 2)
427 {
428 math::sparse_qr<SparseComplexMatrix>
430 if (vector_p)
431 retval = ovl (q.Q (economy), q.R (economy),
432 q.E ().transpose ());
433 else
434 retval = ovl (q.Q (economy), q.R (economy),
435 q.E_MAT ());
436 }
437 else if (nargout > 1)
438 {
439 math::sparse_qr<SparseComplexMatrix>
440 q (arg.sparse_complex_matrix_value (), 0);
441 retval = ovl (q.Q (economy), q.R (economy));
442 }
443 else
444 {
445 math::sparse_qr<SparseComplexMatrix>
446 q (arg.sparse_complex_matrix_value (), 0);
447 retval = ovl (q.R (economy));
448 }
449 }
450 }
451 else
452 {
453 if (have_b && nargout == 1)
454 {
455 octave_idx_type info;
456 if (args(1).issparse () && ! args(1).iscomplex ())
457 retval = ovl (math::sparse_qr<SparseMatrix>::solve
459 (arg.sparse_matrix_value (),
460 args (1).sparse_matrix_value (), info));
461 else if (! args(1).issparse () && args(1).iscomplex ())
462 retval = ovl (math::sparse_qr<SparseMatrix>::solve
464 (arg.sparse_matrix_value (),
465 args (1).complex_matrix_value (), info));
466 else if (! args(1).issparse () && ! args(1).iscomplex ())
467 retval = ovl (math::sparse_qr<SparseMatrix>::solve
469 (arg.sparse_matrix_value (),
470 args (1).matrix_value (), info));
471 else if (args(1).issparse () && args(1).iscomplex ())
472 retval = ovl (math::sparse_qr<SparseMatrix>::solve
474 (arg.sparse_matrix_value (),
475 args(1).sparse_complex_matrix_value (),
476 info));
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 ().transpose ());
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),
505 q.E ().transpose ());
506 else
507 retval = ovl (q.Q (economy), q.R (economy),
508 q.E_MAT ());
509 }
510 else if (nargout > 1)
511 {
512 math::sparse_qr<SparseMatrix>
513 q (arg.sparse_matrix_value (), 0);
514 retval = ovl (q.Q (economy), q.R (economy));
515 }
516 else
517 {
518 math::sparse_qr<SparseMatrix>
519 q (arg.sparse_matrix_value (), 0);
520 retval = ovl (q.R (economy));
521 }
522 }
523 }
524 }
525 else
526 {
527 if (have_b && nargout > 2)
528 error ("qr: too many output arguments when called with A and B");
529
530 if (arg.is_single_type ())
531 {
532 if (arg.isreal ())
533 {
534 math::qr<FloatMatrix>::type type
535 = qr_type<FloatMatrix> (nargout, economy);
536
538
539 switch (nargout)
540 {
541 case 0:
542 case 1:
543 {
544 math::qr<FloatMatrix> fact (m, type);
545 retval = ovl (fact.R ());
546 }
547 break;
548
549 case 2:
550 {
551 math::qr<FloatMatrix> fact (m, type);
552 retval = ovl (fact.Q (), get_qr_r (fact));
553 if (have_b)
554 {
555 if (is_cmplx)
556 retval(0) = fact.Q ().transpose ()
557 * args(1).float_complex_matrix_value ();
558 else
559 retval(0) = fact.Q ().transpose ()
560 * args(1).float_matrix_value ();
561 }
562 }
563 break;
564
565 default:
566 {
567 math::qrp<FloatMatrix> fact (m, type);
568
569 if (vector_p)
570 retval = ovl (fact.Q (), get_qr_r (fact), fact.Pvec ());
571 else
572 retval = ovl (fact.Q (), get_qr_r (fact), fact.P ());
573 }
574 break;
575 }
576 }
577 else if (arg.iscomplex ())
578 {
579 math::qr<FloatComplexMatrix>::type type
580 = qr_type<FloatComplexMatrix> (nargout, economy);
581
583
584 switch (nargout)
585 {
586 case 0:
587 case 1:
588 {
589 math::qr<FloatComplexMatrix> fact (m, type);
590 retval = ovl (fact.R ());
591 }
592 break;
593
594 case 2:
595 {
596 math::qr<FloatComplexMatrix> fact (m, type);
597 retval = ovl (fact.Q (), get_qr_r (fact));
598 if (have_b)
599 retval(0) = conj (fact.Q ().transpose ())
600 * args(1).float_complex_matrix_value ();
601 }
602 break;
603
604 default:
605 {
606 math::qrp<FloatComplexMatrix> fact (m, type);
607 if (vector_p)
608 retval = ovl (fact.Q (), get_qr_r (fact), fact.Pvec ());
609 else
610 retval = ovl (fact.Q (), get_qr_r (fact), fact.P ());
611 }
612 break;
613 }
614 }
615 }
616 else
617 {
618 if (arg.isreal ())
619 {
620 math::qr<Matrix>::type type = 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 (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 (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, "econ");
713%! [qe2, re2] = qr (a, 0);
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, "econ");
724%! [qe2, re2, pe2] = qr (a, 0);
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, "econ");
735%! [qe2, re2] = qr (a, 0);
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, "econ");
746%! [qe2, re2, pe2] = qr (a, 0);
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%!test <*66488>
804%! ## Orientation of 'p' output for dense matrices
805%! [q, r, p] = qr (eye (3));
806%! assert (size (p), [3, 3]);
807%! [q, r, p] = qr (eye (3), 'vector');
808%! assert (size (p), [1, 3]);
809%! [q, r, p] = qr (eye (3), 'econ');
810%! assert (size (p), [3, 3]);
811%! [q, r, p] = qr (eye (3), 0);
812%! assert (size (p), [1, 3]);
813
814%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE") <*66488>
815%! ## Orientation of 'p' output for sparse matrices
816%! [q, r, p] = qr (speye (3));
817%! assert (size (p), [3, 3]);
818%! [q, r, p] = qr (speye (3), 'vector');
819%! assert (size (p), [1, 3]);
820
821%!testif HAVE_SPQR, HAVE_CHOLMOD
822%! [q, r, p] = qr (speye (3), 'econ');
823%! assert (size (p), [3, 3]);
824%! [q, r, p] = qr (speye (3), 0);
825%! assert (size (p), [1, 3]);
826
827## Test input validation
828%!error <Invalid call> qr ()
829%!error <Invalid call> qr (1,2,3,4)
830%!error <too many output arguments> [a,b,c,d] = qr (1)
831%!error <option string must be .*, not "foo"> qr (magic (3), "foo")
832%!error <option string must be .*, not "foo"> qr (magic (3), ones (3, 1), "foo")
833%!error <too many output arguments when called with A and B>
834%! [q, r, p] = qr (ones (3, 2), ones (3, 1));
835%!error <too many output arguments when called with A and B>
836%! [q, r, p] = qr (ones (3, 2), ones (3, 1), 0);
837
838%!function retval = __testqr (q, r, a, p)
839%! tol = 100* eps (class (q));
840%! retval = 0;
841%! if (nargin == 3)
842%! n1 = norm (q*r - a);
843%! n2 = norm (q'*q - eye (columns (q)));
844%! retval = (n1 < tol && n2 < tol);
845%! else
846%! n1 = norm (q'*q - eye (columns (q)));
847%! retval = (n1 < tol);
848%! if (isvector (p))
849%! n2 = norm (q*r - a(:,p));
850%! retval = (retval && n2 < tol);
851%! else
852%! n2 = norm (q*r - a*p);
853%! retval = (retval && n2 < tol);
854%! endif
855%! endif
856%!endfunction
857
858%!test
859%! t = ones (24, 1);
860%! j = 1;
861%!
862%! if (false) # eliminate big matrix tests
863%! a = rand (5000, 20);
864%! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
865%! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
866%! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
867%! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
868%!
869%! a = a+1i*eps;
870%! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
871%! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
872%! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
873%! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
874%! endif
875%!
876%! a = [ ones(1,15); sqrt(eps)*eye(15) ];
877%! [q,r] = qr (a); t(j++) = __testqr (q, r, a);
878%! [q,r] = qr (a'); t(j++) = __testqr (q, r, a');
879%! [q,r,p] = qr (a); t(j++) = __testqr (q, r, a, p);
880%! [q,r,p] = qr (a'); t(j++) = __testqr (q, r, a', p);
881%!
882%! a = a+1i*eps;
883%! [q,r] = qr (a); t(j++) = __testqr (q, r, a);
884%! [q,r] = qr (a'); t(j++) = __testqr (q, r, a');
885%! [q,r,p] = qr (a); t(j++) = __testqr (q, r, a, p);
886%! [q,r,p] = qr (a'); t(j++) = __testqr (q, r, a', p);
887%!
888%! a = [ ones(1,15); sqrt(eps)*eye(15) ];
889%! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
890%! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
891%! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
892%! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
893%!
894%! a = a+1i*eps;
895%! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
896%! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
897%! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
898%! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
899%!
900%! a = [ 611 196 -192 407 -8 -52 -49 29
901%! 196 899 113 -192 -71 -43 -8 -44
902%! -192 113 899 196 61 49 8 52
903%! 407 -192 196 611 8 44 59 -23
904%! -8 -71 61 8 411 -599 208 208
905%! -52 -43 49 44 -599 411 208 208
906%! -49 -8 8 59 208 208 99 -911
907%! 29 -44 52 -23 208 208 -911 99 ];
908%! [q,r] = qr (a);
909%!
910%! assert (all (t));
911%! if (all (t))
912%! assert (norm (q*r - a), 0, 5000*eps);
913%! endif
914
915%!test
916%! a = single ([0, 2, 1; 2, 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, 1; 2, 1, 2]);
926%!
927%! [q, r, p] = qr (a); # FIXME: not giving right dimensions.
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; 2, 1; 1, 2]);
935%!
936%! [q, r] = qr (a);
937%! [qe, re] = qr (a, 0);
938%!
939%! assert (q * r, a, sqrt (eps ("single")));
940%! assert (qe * re, a, sqrt (eps ("single")));
941
942%!test
943%! a = single ([0, 2; 2, 1; 1, 2]);
944%!
945%! [q, r, p] = qr (a);
946%! [qe, re, pe] = qr (a, 0);
947%!
948%! assert (q * r, a * p, sqrt (eps ("single")));
949%! assert (qe * re, a(:, pe), sqrt (eps ("single")));
950
951%!test
952%! a = single ([0, 2, 1; 2, 1, 2; 3, 1, 2]);
953%! b = single ([1, 3, 2; 1, 1, 0; 3, 0, 2]);
954%!
955%! [q, r] = qr (a);
956%! [c, re] = qr (a, b);
957%!
958%! assert (r, re, sqrt (eps ("single")));
959%! assert (q'*b, c, sqrt (eps ("single")));
960
961%!test
962%! a = single ([0, 2, i; 2, 1, 2; 3, 1, 2]);
963%! b = single ([1, 3, 2; 1, i, 0; 3, 0, 2]);
964%!
965%! [q, r] = qr (a);
966%! [c, re] = qr (a, b);
967%!
968%! assert (r, re, sqrt (eps ("single")));
969%! assert (q'*b, c, sqrt (eps ("single")));
970
971%!test
972%! a = single ([0, 2, i; 2, 1, 2; 3, 1, 2]);
973%! b = single ([1, 3, 2; 1, 1, 0; 3, 0, 2]);
974%!
975%! [q, r] = qr (a);
976%! [c, re] = qr (a, b);
977%!
978%! assert (r, re, sqrt (eps));
979%! assert (q'*b, c, sqrt (eps));
980
981%!test
982%! a = single ([0, 2, 1; 2, 1, 2; 3, 1, 2]);
983%! b = single ([1, 3, 2; 1, i, 0; 3, 0, 2]);
984%!
985%! [q, r] = qr (a);
986%! [c, re] = qr (a, b);
987%!
988%! assert (r, re, sqrt (eps ("single")));
989%! assert (q'*b, c, sqrt (eps ("single")));
990
991%!test
992%! t = ones (24, 1);
993%! j = 1;
994%!
995%! if (false) # eliminate big matrix tests
996%! a = rand (5000,20);
997%! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
998%! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
999%! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
1000%! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
1001%!
1002%! a = a+1i* eps ("single");
1003%! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
1004%! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
1005%! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
1006%! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
1007%! endif
1008%!
1009%! a = [ ones(1,15); sqrt(eps("single"))*eye(15) ];
1010%! [q,r] = qr (a); t(j++) = __testqr (q, r, a);
1011%! [q,r] = qr (a'); t(j++) = __testqr (q, r, a');
1012%! [q,r,p] = qr (a); t(j++) = __testqr (q, r, a, p);
1013%! [q,r,p] = qr (a'); t(j++) = __testqr (q, r, a', p);
1014%!
1015%! a = a+1i* eps ("single");
1016%! [q,r] = qr (a); t(j++) = __testqr (q, r, a);
1017%! [q,r] = qr (a'); t(j++) = __testqr (q, r, a');
1018%! [q,r,p] = qr (a); t(j++) = __testqr (q, r, a, p);
1019%! [q,r,p] = qr (a'); t(j++) = __testqr (q, r, a', p);
1020%!
1021%! a = [ ones(1,15); sqrt(eps("single"))*eye(15) ];
1022%! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
1023%! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
1024%! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
1025%! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
1026%!
1027%! a = a+1i* eps ("single");
1028%! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
1029%! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
1030%! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
1031%! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a',p);
1032%!
1033%! a = [ 611 196 -192 407 -8 -52 -49 29
1034%! 196 899 113 -192 -71 -43 -8 -44
1035%! -192 113 899 196 61 49 8 52
1036%! 407 -192 196 611 8 44 59 -23
1037%! -8 -71 61 8 411 -599 208 208
1038%! -52 -43 49 44 -599 411 208 208
1039%! -49 -8 8 59 208 208 99 -911
1040%! 29 -44 52 -23 208 208 -911 99 ];
1041%! [q,r] = qr (a);
1042%!
1043%! assert (all (t));
1044%! if (all (t))
1045%! assert (norm (q*r - a), 0, 5000 * eps ("single"));
1046%! endif
1047
1048## The deactivated tests below can't be tested till rectangular back-subs is
1049## implemented for sparse matrices.
1050
1051%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1052%! n = 20; d = 0.2;
1053%! ## initialize generators to make behavior reproducible
1054%! rand ("state", 42);
1055%! randn ("state", 42);
1056%! a = sprandn (n,n,d) + speye (n,n);
1057%! r = qr (a);
1058%! assert (r'*r, a'*a, 1e-10);
1059
1060%!testif HAVE_COLAMD; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1061%! n = 20; d = 0.2;
1062%! ## initialize generators to make behavior reproducible
1063%! rand ("state", 42);
1064%! randn ("state", 42);
1065%! a = sprandn (n,n,d) + speye (n,n);
1066%! q = symamd (a);
1067%! a = a(q,q);
1068%! r = qr (a);
1069%! assert (r'*r, a'*a, 1e-10);
1070
1071%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1072%! n = 20; d = 0.2;
1073%! ## initialize generators to make behavior reproducible
1074%! rand ("state", 42);
1075%! randn ("state", 42);
1076%! a = sprandn (n,n,d) + speye (n,n);
1077%! [c, r] = qr (a, ones (n,1));
1078%! assert (r\c, full (a)\ones (n,1), 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 = sprandn (n,n,d) + speye (n,n);
1086%! b = randn (n,2);
1087%! [c, r] = qr (a, b);
1088%! assert (r\c, full (a)\b, 10e-10);
1089
1090## Test under-determined systems!!
1091%!#testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1092%! n = 20; d = 0.2;
1093%! ## initialize generators to make behavior reproducible
1094%! rand ("state", 42);
1095%! randn ("state", 42);
1096%! a = sprandn (n,n+1,d) + speye (n,n+1);
1097%! b = randn (n,2);
1098%! [c, r] = qr (a, b);
1099%! assert (r\c, full (a)\b, 10e-10);
1100
1101%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1102%! n = 20; d = 0.2;
1103%! ## initialize generators to make behavior reproducible
1104%! rand ("state", 42);
1105%! randn ("state", 42);
1106%! a = 1i* sprandn (n,n,d) + speye (n,n);
1107%! r = qr (a);
1108%! assert (r'*r, a'*a, 1e-10);
1109
1110%!testif HAVE_COLAMD; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1111%! n = 20; d = 0.2;
1112%! ## initialize generators to make behavior reproducible
1113%! rand ("state", 42);
1114%! randn ("state", 42);
1115%! a = 1i* sprandn (n,n,d) + speye (n,n);
1116%! q = symamd (a);
1117%! a = a(q,q);
1118%! r = qr (a);
1119%! assert (r'*r, a'*a, 1e-10);
1120
1121%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1122%! n = 20; d = 0.2;
1123%! ## initialize generators to make behavior reproducible
1124%! rand ("state", 42);
1125%! randn ("state", 42);
1126%! a = 1i* sprandn (n,n,d) + speye (n,n);
1127%! [c, r] = qr (a, ones (n,1));
1128%! assert (r\c, full (a)\ones (n,1), 10e-10);
1129
1130%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1131%! n = 20; d = 0.2;
1132%! ## initialize generators to make behavior reproducible
1133%! rand ("state", 42);
1134%! randn ("state", 42);
1135%! a = 1i* sprandn (n,n,d) + speye (n,n);
1136%! b = randn (n,2);
1137%! [c, r] = qr (a, b);
1138%! assert (r\c, full (a)\b, 10e-10);
1139
1140## Test under-determined systems!!
1141%!#testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1142%! n = 20; d = 0.2;
1143%! ## initialize generators to make behavior reproducible
1144%! rand ("state", 42);
1145%! randn ("state", 42);
1146%! a = 1i* sprandn (n,n+1,d) + speye (n,n+1);
1147%! b = randn (n, 2);
1148%! [c, r] = qr (a, b);
1149%! assert (r\c, full (a)\b, 10e-10);
1150
1151%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
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] = qr (a, b);
1159%! assert (r\c, full (a)\b, 10e-10);
1160
1161%!testif HAVE_SPQR, HAVE_CHOLMOD
1162%! n = 12; m = 20; d = 0.2;
1163%! ## initialize generators to make behavior reproducible
1164%! rand ("state", 42);
1165%! randn ("state", 42);
1166%! a = sprandn (m, n, d);
1167%! b = sprandn (m, 2, d);
1168%! [c, r] = qr (a, b, 0);
1169%! [c2, r2] = qr (full (a), full (b), 0);
1170%! assert (r\c, r2\c2, 10e-10);
1171
1172%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1173%! n = 12; m = 20; d = 0.2;
1174%! ## initialize generators to make behavior reproducible
1175%! rand ("state", 42);
1176%! randn ("state", 42);
1177%! a = sprandn (m, n, d);
1178%! b = randn (m, 2);
1179%! [c, r, p] = qr (a, b, "matrix");
1180%! x = p * (r\c);
1181%! [c2, r2] = qr (full (a), b);
1182%! x2 = r2 \ c2;
1183%! assert (x, x2, 10e-10);
1184
1185%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1186%! n = 12; m = 20; d = 0.2;
1187%! ## initialize generators to make behavior reproducible
1188%! rand ("state", 42);
1189%! randn ("state", 42);
1190%! a = sprandn (m, n, d);
1191%! [q, r, p] = qr (a, "matrix");
1192%! assert (q * r, a * p, 10e-10);
1193
1194%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1195%! n = 12; m = 20; d = 0.2;
1196%! ## initialize generators to make behavior reproducible
1197%! rand ("state", 42);
1198%! randn ("state", 42);
1199%! a = sprandn (m, n, d);
1200%! b = randn (m, 2);
1201%! x = qr (a, b);
1202%! [c2, r2] = qr (full (a), b);
1203%! assert (x, r2\c2, 10e-10);
1204
1205%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1206%! n = 12; m = 20; d = 0.2;
1207%! ## initialize generators to make behavior reproducible
1208%! rand ("state", 42);
1209%! randn ("state", 42);
1210%! a = sprandn (m, n, d);
1211%! b = i * randn (m, 2);
1212%! x = qr (a, b);
1213%! [c2, r2] = qr (full (a), b);
1214%! assert (x, r2\c2, 10e-10);
1215
1216%!#testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1217%! n = 12; m = 20; d = 0.2;
1218%! ## initialize generators to make behavior reproducible
1219%! rand ("state", 42);
1220%! randn ("state", 42);
1221%! a = sprandn (m, n, d);
1222%! b = i * randn (m, 2);
1223%! [c, r] = qr (a, b);
1224%! [c2, r2] = qr (full (a), b);
1225%! assert (r\c, r2\c2, 10e-10);
1226
1227%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1228%! n = 12; m = 20; d = 0.2;
1229%! ## initialize generators to make behavior reproducible
1230%! rand ("state", 42);
1231%! randn ("state", 42);
1232%! a = sprandn (m, n, d);
1233%! b = i * randn (m, 2);
1234%! [c, r, p] = qr (a, b, "matrix");
1235%! x = p * (r\c);
1236%! [c2, r2] = qr (full (a), b);
1237%! x2 = r2 \ c2;
1238%! assert (x, x2, 10e-10);
1239
1240%!testif HAVE_SPQR, HAVE_CHOLMOD
1241%! n = 12; m = 20; d = 0.2;
1242%! ## initialize generators to make behavior reproducible
1243%! rand ("state", 42);
1244%! randn ("state", 42);
1245%! a = i * sprandn (m, n, d);
1246%! b = sprandn (m, 2, d);
1247%! [c, r] = qr (a, b, 0);
1248%! [c2, r2] = qr (full (a), full (b), 0);
1249%! assert (r\c, r2\c2, 10e-10);
1250
1251%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1252%! n = 12; m = 20; d = 0.2;
1253%! ## initialize generators to make behavior reproducible
1254%! rand ("state", 42);
1255%! randn ("state", 42);
1256%! a = i * sprandn (m, n, d);
1257%! b = randn (m, 2);
1258%! [c, r, p] = qr (a, b, "matrix");
1259%! x = p * (r\c);
1260%! [c2, r2] = qr (full (a), b);
1261%! x2 = r2 \ c2;
1262%! assert(x, x2, 10e-10);
1263
1264%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1265%! n = 12; m = 20; d = 0.2;
1266%! ## initialize generators to make behavior reproducible
1267%! rand ("state", 42);
1268%! randn ("state", 42);
1269%! a = i * sprandn (m, n, d);
1270%! [q, r, p] = qr (a, "matrix");
1271%! assert(q * r, a * p, 10e-10);
1272
1273%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1274%! n = 12; m = 20; d = 0.2;
1275%! ## initialize generators to make behavior reproducible
1276%! rand ("state", 42);
1277%! randn ("state", 42);
1278%! a = i * sprandn (m, n, d);
1279%! b = randn (m, 2);
1280%! x = qr (a, b);
1281%! [c2, r2] = qr (full (a), b);
1282%! assert (x, r2\c2, 10e-10);
1283
1284%!testif HAVE_SPQR, HAVE_CHOLMOD
1285%! a = sparse (5, 6);
1286%! a(3,1) = 0.8;
1287%! a(2,2) = 1.4;
1288%! a(1,6) = -0.5;
1289%! r = qr (a);
1290%! assert (r'*r, a'*a, 10e-10);
1291*/
1292
1293static
1294bool
1295check_qr_dims (const octave_value& q, const octave_value& r,
1296 bool allow_ecf = false)
1297{
1298 octave_idx_type m = q.rows ();
1299 octave_idx_type k = r.rows ();
1300 octave_idx_type n = r.columns ();
1301 return ((q.ndims () == 2 && r.ndims () == 2 && k == q.columns ())
1302 && (m == k || (allow_ecf && k == n && k < m)));
1303}
1304
1305static
1306bool
1307check_index (const octave_value& i, bool vector_allowed = false)
1308{
1309 return ((i.isreal () || i.isinteger ())
1310 && (i.is_scalar_type () || vector_allowed));
1311}
1312
1313DEFUN (qrupdate, args, ,
1314 doc: /* -*- texinfo -*-
1315@deftypefn {} {[@var{Q1}, @var{R1}] =} qrupdate (@var{Q}, @var{R}, @var{u}, @var{v})
1316Update a QR factorization given update vectors or matrices.
1317
1318Given a QR@tie{}factorization of a real or complex matrix
1319@w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and
1320@var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization of
1321@w{@var{A} + @var{u}*@var{v}'}, where @var{u} and @var{v} are column vectors
1322(rank-1 update) or matrices with equal number of columns
1323(rank-k update). Notice that the latter case is done as a sequence of
1324rank-1 updates; thus, for k large enough, it will be both faster and more
1325accurate to recompute the factorization from scratch.
1326
1327The QR@tie{}factorization supplied may be either full (Q is square) or
1328economized (R is square).
1329
1330@seealso{qr, qrinsert, qrdelete, qrshift}
1331@end deftypefn */)
1332{
1333 octave_value_list retval;
1334
1335 if (args.length () != 4)
1336 print_usage ();
1337
1338 octave_value argq = args(0);
1339 octave_value argr = args(1);
1340 octave_value argu = args(2);
1341 octave_value argv = args(3);
1342
1343 if (! argq.isnumeric () || ! argr.isnumeric ()
1344 || ! argu.isnumeric () || ! argv.isnumeric ())
1345 print_usage ();
1346
1347 if (! check_qr_dims (argq, argr, true))
1348 error ("qrupdate: Q and R dimensions don't match");
1349
1350 if (argq.isreal () && argr.isreal () && argu.isreal ()
1351 && argv.isreal ())
1352 {
1353 // all real case
1354 if (argq.is_single_type () || argr.is_single_type ()
1355 || argu.is_single_type () || argv.is_single_type ())
1356 {
1358 FloatMatrix R = argr.float_matrix_value ();
1359 FloatMatrix u = argu.float_matrix_value ();
1360 FloatMatrix v = argv.float_matrix_value ();
1361
1362 math::qr<FloatMatrix> fact (Q, R);
1363 fact.update (u, v);
1364
1365 retval = ovl (fact.Q (), get_qr_r (fact));
1366 }
1367 else
1368 {
1369 Matrix Q = argq.matrix_value ();
1370 Matrix R = argr.matrix_value ();
1371 Matrix u = argu.matrix_value ();
1372 Matrix v = argv.matrix_value ();
1373
1374 math::qr<Matrix> fact (Q, R);
1375 fact.update (u, v);
1376
1377 retval = ovl (fact.Q (), get_qr_r (fact));
1378 }
1379 }
1380 else
1381 {
1382 // complex case
1383 if (argq.is_single_type () || argr.is_single_type ()
1384 || argu.is_single_type () || argv.is_single_type ())
1385 {
1390
1391 math::qr<FloatComplexMatrix> fact (Q, R);
1392 fact.update (u, v);
1393
1394 retval = ovl (fact.Q (), get_qr_r (fact));
1395 }
1396 else
1397 {
1402
1403 math::qr<ComplexMatrix> fact (Q, R);
1404 fact.update (u, v);
1405
1406 retval = ovl (fact.Q (), get_qr_r (fact));
1407 }
1408 }
1409
1410 return retval;
1411}
1412
1413/*
1414%!shared A, u, v, Ac, uc, vc
1415%! A = [0.091364 0.613038 0.999083;
1416%! 0.594638 0.425302 0.603537;
1417%! 0.383594 0.291238 0.085574;
1418%! 0.265712 0.268003 0.238409;
1419%! 0.669966 0.743851 0.445057 ];
1420%!
1421%! u = [0.85082;
1422%! 0.76426;
1423%! 0.42883;
1424%! 0.53010;
1425%! 0.80683 ];
1426%!
1427%! v = [0.98810;
1428%! 0.24295;
1429%! 0.43167 ];
1430%!
1431%! Ac = [0.620405 + 0.956953i 0.480013 + 0.048806i 0.402627 + 0.338171i;
1432%! 0.589077 + 0.658457i 0.013205 + 0.279323i 0.229284 + 0.721929i;
1433%! 0.092758 + 0.345687i 0.928679 + 0.241052i 0.764536 + 0.832406i;
1434%! 0.912098 + 0.721024i 0.049018 + 0.269452i 0.730029 + 0.796517i;
1435%! 0.112849 + 0.603871i 0.486352 + 0.142337i 0.355646 + 0.151496i ];
1436%!
1437%! uc = [0.20351 + 0.05401i;
1438%! 0.13141 + 0.43708i;
1439%! 0.29808 + 0.08789i;
1440%! 0.69821 + 0.38844i;
1441%! 0.74871 + 0.25821i ];
1442%!
1443%! vc = [0.85839 + 0.29468i;
1444%! 0.20820 + 0.93090i;
1445%! 0.86184 + 0.34689i ];
1446%!
1447
1448%!test
1449%! [Q,R] = qr (A);
1450%! [Q,R] = qrupdate (Q, R, u, v);
1451%! assert (norm (vec (Q'*Q - eye (5)), Inf), 0, 1e1*eps);
1452%! assert (norm (vec (triu (R)-R), Inf), 0);
1453%! assert (norm (vec (Q*R - A - u*v'), Inf), 0, norm (A)*1e1*eps);
1454%!
1455%!test
1456%! [Q,R] = qr (Ac);
1457%! [Q,R] = qrupdate (Q, R, uc, vc);
1458%! assert (norm (vec (Q'*Q - eye (5)), Inf), 0, 1e1*eps);
1459%! assert (norm (vec (triu (R)-R), Inf), 0);
1460%! assert (norm (vec (Q*R - Ac - uc*vc'), Inf), 0, norm (Ac)*1e1*eps);
1461
1462%!test
1463%! [Q,R] = qr (single (A));
1464%! [Q,R] = qrupdate (Q, R, single (u), single (v));
1465%! assert (norm (vec (Q'*Q - eye (5,"single")), Inf), single (0), ...
1466%! 1e1 * eps ("single"));
1467%! assert (norm (vec (triu (R)-R), Inf), single (0));
1468%! assert (norm (vec (Q*R - single (A) - single (u)* single (v)'), Inf), ...
1469%! single (0), norm (single (A))*1e1 * eps ("single"));
1470%!
1471%!test
1472%! [Q,R] = qr (single (Ac));
1473%! [Q,R] = qrupdate (Q, R, single (uc), single (vc));
1474%! assert (norm (vec (Q'*Q - eye (5,"single")), Inf), single (0), ...
1475%! 1e1 * eps ("single"));
1476%! assert (norm (vec (triu (R)-R), Inf), single (0));
1477%! assert (norm (vec (Q*R - single (Ac) - single (uc)* single (vc)'), Inf), ...
1478%! single (0), norm (single (Ac))*1e1 * eps ("single"));
1479*/
1480
1481DEFUN (qrinsert, args, ,
1482 doc: /* -*- texinfo -*-
1483@deftypefn {} {[@var{Q1}, @var{R1}] =} qrinsert (@var{Q}, @var{R}, @var{j}, @var{x}, @var{orient})
1484Update a QR factorization given a row or column to insert in the original
1485factored matrix.
1486
1487
1488Given a QR@tie{}factorization of a real or complex matrix
1489@w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and
1490@var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization of
1491@w{[A(:,1:j-1) x A(:,j:n)]}, where @var{u} is a column vector to be inserted
1492into @var{A} (if @var{orient} is @qcode{"col"}), or the
1493QR@tie{}factorization of @w{[A(1:j-1,:);x;A(:,j:n)]}, where @var{x} is a row
1494vector to be inserted into @var{A} (if @var{orient} is @qcode{"row"}).
1495
1496The default value of @var{orient} is @qcode{"col"}. If @var{orient} is
1497@qcode{"col"}, @var{u} may be a matrix and @var{j} an index vector
1498resulting in the QR@tie{}factorization of a matrix @var{B} such that
1499@w{B(:,@var{j})}@ gives @var{u} and @w{B(:,@var{j}) = []}@ gives @var{A}.
1500Notice that the latter case is done as a sequence of k insertions;
1501thus, for k large enough, it will be both faster and more accurate to
1502recompute the factorization from scratch.
1503
1504If @var{orient} is @qcode{"col"}, the QR@tie{}factorization supplied may
1505be either full (Q is square) or economized (R is square).
1506
1507If @var{orient} is @qcode{"row"}, full factorization is needed.
1508@seealso{qr, qrupdate, qrdelete, qrshift}
1509@end deftypefn */)
1510{
1511 int nargin = args.length ();
1512
1513 if (nargin < 4 || nargin > 5)
1514 print_usage ();
1515
1516 octave_value argq = args(0);
1517 octave_value argr = args(1);
1518 octave_value argj = args(2);
1519 octave_value argx = args(3);
1520
1521 if (! argq.isnumeric () || ! argr.isnumeric ()
1522 || ! argx.isnumeric ()
1523 || (nargin > 4 && ! args(4).is_string ()))
1524 print_usage ();
1525
1526 std::string orient = (nargin < 5) ? "col" : args(4).string_value ();
1527 bool col = (orient == "col");
1528
1529 if (! col && orient != "row")
1530 error (R"(qrinsert: ORIENT must be "col" or "row")");
1531
1532 if (! check_qr_dims (argq, argr, col) || (! col && argx.rows () != 1))
1533 error ("qrinsert: dimension mismatch");
1534
1535 if (! check_index (argj, col))
1536 error ("qrinsert: invalid index J");
1537
1538 octave_value_list retval;
1539
1541
1542 octave_idx_type one = 1;
1543
1544 if (argq.isreal () && argr.isreal () && argx.isreal ())
1545 {
1546 // real case
1547 if (argq.is_single_type () || argr.is_single_type ()
1548 || argx.is_single_type ())
1549 {
1551 FloatMatrix R = argr.float_matrix_value ();
1553
1554 math::qr<FloatMatrix> fact (Q, R);
1555
1556 if (col)
1557 fact.insert_col (x, j-one);
1558 else
1559 fact.insert_row (x.row (0), j(0)-one);
1560
1561 retval = ovl (fact.Q (), get_qr_r (fact));
1562 }
1563 else
1564 {
1565 Matrix Q = argq.matrix_value ();
1566 Matrix R = argr.matrix_value ();
1567 Matrix x = argx.matrix_value ();
1568
1569 math::qr<Matrix> fact (Q, R);
1570
1571 if (col)
1572 fact.insert_col (x, j-one);
1573 else
1574 fact.insert_row (x.row (0), j(0)-one);
1575
1576 retval = ovl (fact.Q (), get_qr_r (fact));
1577 }
1578 }
1579 else
1580 {
1581 // complex case
1582 if (argq.is_single_type () || argr.is_single_type ()
1583 || argx.is_single_type ())
1584 {
1588
1589 math::qr<FloatComplexMatrix> fact (Q, R);
1590
1591 if (col)
1592 fact.insert_col (x, j-one);
1593 else
1594 fact.insert_row (x.row (0), j(0)-one);
1595
1596 retval = ovl (fact.Q (), get_qr_r (fact));
1597 }
1598 else
1599 {
1603
1604 math::qr<ComplexMatrix> fact (Q, R);
1605
1606 if (col)
1607 fact.insert_col (x, j-one);
1608 else
1609 fact.insert_row (x.row (0), j(0)-one);
1610
1611 retval = ovl (fact.Q (), get_qr_r (fact));
1612 }
1613 }
1614
1615 return retval;
1616}
1617
1618/*
1619%!test
1620%! [Q,R] = qr (A);
1621%! [Q,R] = qrinsert (Q, R, 3, u);
1622%! assert (norm (vec (Q'*Q - eye (5)), Inf), 0, 1e1*eps);
1623%! assert (norm (vec (triu (R) - R), Inf), 0);
1624%! assert (norm (vec (Q*R - [A(:,1:2) u A(:,3)]), Inf), 0, norm (A)*1e1*eps);
1625%!test
1626%! [Q,R] = qr (Ac);
1627%! [Q,R] = qrinsert (Q, R, 3, uc);
1628%! assert (norm (vec (Q'*Q - eye (5)), Inf), 0, 1e1*eps);
1629%! assert (norm (vec (triu (R) - R), Inf), 0);
1630%! assert (norm (vec (Q*R - [Ac(:,1:2) uc Ac(:,3)]), Inf), 0, ...
1631%! norm (Ac) * 1e1 * eps);
1632%!test
1633%! x = [0.85082 0.76426 0.42883 ];
1634%!
1635%! [Q,R] = qr (A);
1636%! [Q,R] = qrinsert (Q, R, 3, x, "row");
1637%! assert (norm (vec (Q'*Q - eye (6)), Inf), 0, 1e1*eps);
1638%! assert (norm (vec (triu (R) - R), Inf), 0);
1639%! assert (norm (vec (Q*R - [A(1:2,:);x;A(3:5,:)]), Inf), 0, norm (A)*1e1*eps);
1640%!test
1641%! x = [0.20351 + 0.05401i 0.13141 + 0.43708i 0.29808 + 0.08789i ];
1642%!
1643%! [Q,R] = qr (Ac);
1644%! [Q,R] = qrinsert (Q, R, 3, x, "row");
1645%! assert (norm (vec (Q'*Q - eye (6)), Inf), 0, 1e1*eps);
1646%! assert (norm (vec (triu (R) - R), Inf), 0);
1647%! assert (norm (vec (Q*R - [Ac(1:2,:);x;Ac(3:5,:)]), Inf), 0, ...
1648%! norm (Ac) * 1e1 * eps);
1649
1650%!test
1651%! [Q,R] = qr (single (A));
1652%! [Q,R] = qrinsert (Q, R, 3, single (u));
1653%! assert (norm (vec (Q'*Q - eye (5,"single")), Inf), single (0), ...
1654%! 1e1 * eps ("single"));
1655%! assert (norm (vec (triu (R) - R), Inf), single (0));
1656%! assert (norm (vec (Q*R - single ([A(:,1:2) u A(:,3)])), Inf), ...
1657%! single (0), norm (single (A))*1e1 * eps ("single"));
1658%!test
1659%! [Q,R] = qr (single (Ac));
1660%! [Q,R] = qrinsert (Q, R, 3, single (uc));
1661%! assert (norm (vec (Q'*Q - eye (5,"single")), Inf), single (0), ...
1662%! 1e1 * eps ("single"));
1663%! assert (norm (vec (triu (R) - R), Inf), single (0));
1664%! assert (norm (vec (Q*R - single ([Ac(:,1:2) uc Ac(:,3)])), Inf), ...
1665%! single (0), norm (single (Ac))*1e1 * eps ("single"));
1666%!test
1667%! x = single ([0.85082 0.76426 0.42883 ]);
1668%!
1669%! [Q,R] = qr (single (A));
1670%! [Q,R] = qrinsert (Q, R, 3, x, "row");
1671%! assert (norm (vec (Q'*Q - eye (6,"single")), Inf), single (0), ...
1672%! 1e1 * eps ("single"));
1673%! assert (norm (vec (triu (R) - R), Inf), single (0));
1674%! assert (norm (vec (Q*R - single ([A(1:2,:);x;A(3:5,:)])), Inf), ...
1675%! single (0), norm (single (A))*1e1 * eps ("single"));
1676%!test
1677%! x = single ([0.20351 + 0.05401i 0.13141 + 0.43708i 0.29808 + 0.08789i ]);
1678%!
1679%! [Q,R] = qr (single (Ac));
1680%! [Q,R] = qrinsert (Q, R, 3, x, "row");
1681%! assert (norm (vec (Q'*Q - eye (6,"single")), Inf), single (0), ...
1682%! 1e1 * eps ("single"));
1683%! assert (norm (vec (triu (R) - R), Inf), single (0));
1684%! assert (norm (vec (Q*R - single ([Ac(1:2,:);x;Ac(3:5,:)])), Inf), ...
1685%! single (0), norm (single (Ac))*1e1 * eps ("single"));
1686*/
1687
1688DEFUN (qrdelete, args, ,
1689 doc: /* -*- texinfo -*-
1690@deftypefn {} {[@var{Q1}, @var{R1}] =} qrdelete (@var{Q}, @var{R}, @var{j}, @var{orient})
1691Update a QR factorization given a row or column to delete from the original
1692factored matrix.
1693
1694Given a QR@tie{}factorization of a real or complex matrix
1695@w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and
1696@var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization of
1697@w{[A(:,1:j-1), U, A(:,j:n)]},
1698where @var{u} is a column vector to be inserted into @var{A}
1699(if @var{orient} is @qcode{"col"}),
1700or the QR@tie{}factorization of @w{[A(1:j-1,:);X;A(:,j:n)]},
1701where @var{x} is a row @var{orient} is @qcode{"row"}).
1702The default value of @var{orient} is @qcode{"col"}.
1703
1704If @var{orient} is @qcode{"col"}, @var{j} may be an index vector
1705resulting in the QR@tie{}factorization of a matrix @var{B} such that
1706@w{A(:,@var{j}) = []}@ gives @var{B}. Notice that the latter case is done as
1707a sequence of k deletions; thus, for k large enough, it will be both faster
1708and more accurate to recompute the factorization from scratch.
1709
1710If @var{orient} is @qcode{"col"}, the QR@tie{}factorization supplied may
1711be either full (Q is square) or economized (R is square).
1712
1713If @var{orient} is @qcode{"row"}, full factorization is needed.
1714@seealso{qr, qrupdate, qrinsert, qrshift}
1715@end deftypefn */)
1716{
1717 int nargin = args.length ();
1718
1719 if (nargin < 3 || nargin > 4)
1720 print_usage ();
1721
1722 octave_value argq = args(0);
1723 octave_value argr = args(1);
1724 octave_value argj = args(2);
1725
1726 if (! argq.isnumeric () || ! argr.isnumeric ()
1727 || (nargin > 3 && ! args(3).is_string ()))
1728 print_usage ();
1729
1730 std::string orient = (nargin < 4) ? "col" : args(3).string_value ();
1731 bool col = orient == "col";
1732
1733 if (! col && orient != "row")
1734 error (R"(qrdelete: ORIENT must be "col" or "row")");
1735
1736 if (! check_qr_dims (argq, argr, col))
1737 error ("qrdelete: dimension mismatch");
1738
1740 if (! check_index (argj, col))
1741 error ("qrdelete: invalid index J");
1742
1743 octave_value_list retval;
1744
1745 octave_idx_type one = 1;
1746
1747 if (argq.isreal () && argr.isreal ())
1748 {
1749 // real case
1750 if (argq.is_single_type () || argr.is_single_type ())
1751 {
1753 FloatMatrix R = argr.float_matrix_value ();
1754
1755 math::qr<FloatMatrix> fact (Q, R);
1756
1757 if (col)
1758 fact.delete_col (j-one);
1759 else
1760 fact.delete_row (j(0)-one);
1761
1762 retval = ovl (fact.Q (), get_qr_r (fact));
1763 }
1764 else
1765 {
1766 Matrix Q = argq.matrix_value ();
1767 Matrix R = argr.matrix_value ();
1768
1769 math::qr<Matrix> fact (Q, R);
1770
1771 if (col)
1772 fact.delete_col (j-one);
1773 else
1774 fact.delete_row (j(0)-one);
1775
1776 retval = ovl (fact.Q (), get_qr_r (fact));
1777 }
1778 }
1779 else
1780 {
1781 // complex case
1782 if (argq.is_single_type () || argr.is_single_type ())
1783 {
1786
1787 math::qr<FloatComplexMatrix> fact (Q, R);
1788
1789 if (col)
1790 fact.delete_col (j-one);
1791 else
1792 fact.delete_row (j(0)-one);
1793
1794 retval = ovl (fact.Q (), get_qr_r (fact));
1795 }
1796 else
1797 {
1800
1801 math::qr<ComplexMatrix> fact (Q, R);
1802
1803 if (col)
1804 fact.delete_col (j-one);
1805 else
1806 fact.delete_row (j(0)-one);
1807
1808 retval = ovl (fact.Q (), get_qr_r (fact));
1809 }
1810 }
1811
1812 return retval;
1813}
1814
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);
1825%! assert (norm (vec (Q'*Q - eye (5)), Inf), 0, 16*eps);
1826%! assert (norm (vec (triu (R) - R), Inf), 0);
1827%! assert (norm (vec (Q*R - [AA(:,1:2) AA(:,4)]), 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);
1838%! assert (norm (vec (Q'*Q - eye (5)), Inf) < 16*eps);
1839%! assert (norm (vec (triu (R) - R), Inf) == 0);
1840%! assert (norm (vec (Q*R - [AA(:,1:2) AA(:,4)]), Inf) < norm (AA)*1e1*eps);
1841%!
1842%!test
1843%! AA = [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, "row");
1851%! assert (norm (vec (Q'*Q - eye (4)), Inf), 0, 1e1*eps);
1852%! assert (norm (vec (triu (R) - R), Inf), 0);
1853%! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf), 0, norm (AA)*1e1*eps);
1854%!
1855%!test
1856%! AA = [0.364554 + 0.993117i 0.669818 + 0.510234i 0.426568 + 0.041337i 0.847051 + 0.233291i;
1857%! 0.049600 + 0.242783i 0.448946 + 0.484022i 0.141155 + 0.074420i 0.446746 + 0.392706i;
1858%! 0.581922 + 0.657416i 0.581460 + 0.030016i 0.219909 + 0.447288i 0.201144 + 0.069132i;
1859%! 0.694986 + 0.000571i 0.682327 + 0.841712i 0.807537 + 0.166086i 0.192767 + 0.358098i;
1860%! 0.945002 + 0.066788i 0.350492 + 0.642638i 0.579629 + 0.048102i 0.600170 + 0.636938i ] * I;
1861%!
1862%! [Q,R] = qr (AA);
1863%! [Q,R] = qrdelete (Q, R, 3, "row");
1864%! assert (norm (vec (Q'*Q - eye (4)), Inf), 0, 1e1*eps);
1865%! assert (norm (vec (triu (R) - R), Inf), 0);
1866%! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf), 0, norm (AA)*1e1*eps);
1867
1868%!test
1869%! AA = single ([0.091364 0.613038 0.027504 0.999083;
1870%! 0.594638 0.425302 0.562834 0.603537;
1871%! 0.383594 0.291238 0.742073 0.085574;
1872%! 0.265712 0.268003 0.783553 0.238409;
1873%! 0.669966 0.743851 0.457255 0.445057 ]);
1874%!
1875%! [Q,R] = qr (AA);
1876%! [Q,R] = qrdelete (Q, R, 3);
1877%! assert (norm (vec (Q'*Q - eye (5,"single")), Inf), single (0), ...
1878%! 1e1 * eps ("single"));
1879%! assert (norm (vec (triu (R) - R), Inf), single (0));
1880%! assert (norm (vec (Q*R - [AA(:,1:2) AA(:,4)]), Inf), single (0), ...
1881%! norm (AA)*1e1 * eps ("single"));
1882%!
1883%!test
1884%! AA = single ([0.364554 + 0.993117i 0.669818 + 0.510234i 0.426568 + 0.041337i 0.847051 + 0.233291i;
1885%! 0.049600 + 0.242783i 0.448946 + 0.484022i 0.141155 + 0.074420i 0.446746 + 0.392706i;
1886%! 0.581922 + 0.657416i 0.581460 + 0.030016i 0.219909 + 0.447288i 0.201144 + 0.069132i;
1887%! 0.694986 + 0.000571i 0.682327 + 0.841712i 0.807537 + 0.166086i 0.192767 + 0.358098i;
1888%! 0.945002 + 0.066788i 0.350492 + 0.642638i 0.579629 + 0.048102i 0.600170 + 0.636938i ]) * I;
1889%!
1890%! [Q,R] = qr (AA);
1891%! [Q,R] = qrdelete (Q, R, 3);
1892%! assert (norm (vec (Q'*Q - eye (5,"single")), Inf), single (0), ...
1893%! 1e1 * eps ("single"));
1894%! assert (norm (vec (triu (R) - R), Inf), single (0));
1895%! assert (norm (vec (Q*R - [AA(:,1:2) AA(:,4)]), Inf), single (0), ...
1896%! norm (AA)*1e1 * eps ("single"));
1897
1898%!test
1899%! AA = single ([0.091364 0.613038 0.027504 0.999083;
1900%! 0.594638 0.425302 0.562834 0.603537;
1901%! 0.383594 0.291238 0.742073 0.085574;
1902%! 0.265712 0.268003 0.783553 0.238409;
1903%! 0.669966 0.743851 0.457255 0.445057 ]);
1904%!
1905%! [Q,R] = qr (AA);
1906%! [Q,R] = qrdelete (Q, R, 3, "row");
1907%! assert (norm (vec (Q'*Q - eye (4,"single")), Inf), single (0), ...
1908%! 1.5e1 * eps ("single"));
1909%! assert (norm (vec (triu (R) - R), Inf), single (0));
1910%! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf), single (0), ...
1911%! norm (AA)*1e1 * eps ("single"));
1912%!testif HAVE_QRUPDATE
1913%! ## Same test as above but with more precicision
1914%! AA = single ([0.091364 0.613038 0.027504 0.999083;
1915%! 0.594638 0.425302 0.562834 0.603537;
1916%! 0.383594 0.291238 0.742073 0.085574;
1917%! 0.265712 0.268003 0.783553 0.238409;
1918%! 0.669966 0.743851 0.457255 0.445057 ]);
1919%!
1920%! [Q,R] = qr (AA);
1921%! [Q,R] = qrdelete (Q, R, 3, "row");
1922%! assert (norm (vec (Q'*Q - eye (4,"single")), Inf), single (0), ...
1923%! 1e1* eps ("single"));
1924%! assert (norm (vec (triu (R) - R), Inf), single (0));
1925%! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf), single (0), ...
1926%! norm (AA)*1e1 * eps ("single"));
1927%!
1928%!test
1929%! AA = single ([0.364554 + 0.993117i 0.669818 + 0.510234i 0.426568 + 0.041337i 0.847051 + 0.233291i;
1930%! 0.049600 + 0.242783i 0.448946 + 0.484022i 0.141155 + 0.074420i 0.446746 + 0.392706i;
1931%! 0.581922 + 0.657416i 0.581460 + 0.030016i 0.219909 + 0.447288i 0.201144 + 0.069132i;
1932%! 0.694986 + 0.000571i 0.682327 + 0.841712i 0.807537 + 0.166086i 0.192767 + 0.358098i;
1933%! 0.945002 + 0.066788i 0.350492 + 0.642638i 0.579629 + 0.048102i 0.600170 + 0.636938i ]) * I;
1934%!
1935%! [Q,R] = qr (AA);
1936%! [Q,R] = qrdelete (Q, R, 3, "row");
1937%! assert (norm (vec (Q'*Q - eye (4,"single")), Inf), single (0), ...
1938%! 1e1 * eps ("single"));
1939%! assert (norm (vec (triu (R) - R), Inf), single (0));
1940%! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf), single (0), ...
1941%! norm (AA)*1e1 * eps ("single"));
1942*/
1943
1944DEFUN (qrshift, args, ,
1945 doc: /* -*- texinfo -*-
1946@deftypefn {} {[@var{Q1}, @var{R1}] =} qrshift (@var{Q}, @var{R}, @var{i}, @var{j})
1947Update a QR factorization given a range of columns to shift in the original
1948factored matrix.
1949
1950Given a QR@tie{}factorization of a real or complex matrix
1951@w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and
1952@var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization
1953of @w{@var{A}(:,p)}, where @w{p} is the permutation @*
1954@code{p = [1:i-1, circshift(i:j, 1), j+1:n]} if @w{@var{i} < @var{j}} @*
1955 or @*
1956@code{p = [1:j-1, circshift(j:i,-1), i+1:n]} if @w{@var{j} < @var{i}}. @*
1957
1958@seealso{qr, qrupdate, qrinsert, qrdelete}
1959@end deftypefn */)
1960{
1961 if (args.length () != 4)
1962 print_usage ();
1963
1964 octave_value argq = args(0);
1965 octave_value argr = args(1);
1966 octave_value argi = args(2);
1967 octave_value argj = args(3);
1968
1969 if (! argq.isnumeric () || ! argr.isnumeric ())
1970 print_usage ();
1971
1972 if (! check_qr_dims (argq, argr, true))
1973 error ("qrshift: dimensions mismatch");
1974
1975 octave_idx_type i = argi.idx_type_value ();
1976 octave_idx_type j = argj.idx_type_value ();
1977
1978 if (! check_index (argi) || ! check_index (argj))
1979 error ("qrshift: invalid index I or J");
1980
1981 octave_value_list retval;
1982
1983 if (argq.isreal () && argr.isreal ())
1984 {
1985 // all real case
1986 if (argq.is_single_type ()
1987 && argr.is_single_type ())
1988 {
1990 FloatMatrix R = argr.float_matrix_value ();
1991
1992 math::qr<FloatMatrix> fact (Q, R);
1993 fact.shift_cols (i-1, j-1);
1994
1995 retval = ovl (fact.Q (), get_qr_r (fact));
1996 }
1997 else
1998 {
1999 Matrix Q = argq.matrix_value ();
2000 Matrix R = argr.matrix_value ();
2001
2002 math::qr<Matrix> fact (Q, R);
2003 fact.shift_cols (i-1, j-1);
2004
2005 retval = ovl (fact.Q (), get_qr_r (fact));
2006 }
2007 }
2008 else
2009 {
2010 // complex case
2011 if (argq.is_single_type ()
2012 && argr.is_single_type ())
2013 {
2016
2017 math::qr<FloatComplexMatrix> fact (Q, R);
2018 fact.shift_cols (i-1, j-1);
2019
2020 retval = ovl (fact.Q (), get_qr_r (fact));
2021 }
2022 else
2023 {
2026
2027 math::qr<ComplexMatrix> fact (Q, R);
2028 fact.shift_cols (i-1, j-1);
2029
2030 retval = ovl (fact.Q (), get_qr_r (fact));
2031 }
2032 }
2033
2034 return retval;
2035}
2036
2037/*
2038%!test
2039%! AA = A.';
2040%! i = 2; j = 4; p = [1:i-1, circshift(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)), 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%! j = 2; i = 4; p = [1:j-1, circshift(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)), Inf), 0, 1e1*eps);
2053%! assert (norm (vec (triu (R) - R), Inf), 0);
2054%! assert (norm (vec (Q*R - AA(:,p)), Inf), 0, norm (AA)*1e1*eps);
2055%!
2056%!test
2057%! AA = Ac.';
2058%! i = 2; j = 4; p = [1:i-1, circshift(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)), Inf), 0, 1e1*eps);
2063%! assert (norm (vec (triu (R) - R), Inf), 0);
2064%! assert (norm (vec (Q*R - AA(:,p)), Inf), 0, norm (AA)*1e1*eps);
2065%!
2066%! j = 2; i = 4; p = [1:j-1, circshift(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)), Inf), 0, 1e1*eps);
2071%! assert (norm (vec (triu (R) - R), Inf), 0);
2072%! assert (norm (vec (Q*R - AA(:,p)), Inf), 0, norm (AA)*1e1*eps);
2073
2074%!test
2075%! AA = single (A).';
2076%! i = 2; j = 4; p = [1:i-1, circshift(i:j,-1), j+1:5];
2077%!
2078%! [Q,R] = qr (AA);
2079%! [Q,R] = qrshift (Q, R, i, j);
2080%! assert (norm (vec (Q'*Q - eye (3,"single")), Inf), single (0), ...
2081%! 1e1 * eps ("single"));
2082%! assert (norm (vec (triu (R) - R), Inf), single (0));
2083%! assert (norm (vec (Q*R - AA(:,p)), Inf), single (0), ...
2084%! norm (AA)*1e1 * eps ("single"));
2085%!
2086%! j = 2; i = 4; p = [1:j-1, circshift(j:i,+1), i+1:5];
2087%!
2088%! [Q,R] = qr (AA);
2089%! [Q,R] = qrshift (Q, R, i, j);
2090%! assert (norm (vec (Q'*Q - eye (3,"single")), Inf), single (0), ...
2091%! 1e1 * eps ("single"));
2092%! assert (norm (vec (triu (R) - R), Inf), single (0));
2093%! assert (norm (vec (Q*R - AA(:,p)), Inf), single (0), ...
2094%! norm (AA)*1e1 * eps ("single"));
2095%!
2096%!test
2097%! AA = single (Ac).';
2098%! i = 2; j = 4; p = [1:i-1, circshift(i:j,-1), j+1:5];
2099%!
2100%! [Q,R] = qr (AA);
2101%! [Q,R] = qrshift (Q, R, i, j);
2102%! assert (norm (vec (Q'*Q - eye (3,"single")), Inf), single (0), ...
2103%! 1e1 * eps ("single"));
2104%! assert (norm (vec (triu (R) - R), Inf), single (0));
2105%! assert (norm (vec (Q*R - AA(:,p)), Inf), single (0), ...
2106%! norm (AA)*1e1 * eps ("single"));
2107%!
2108%! j = 2; i = 4; p = [1:j-1, circshift(j:i,+1), i+1:5];
2109%!
2110%! [Q,R] = qr (AA);
2111%! [Q,R] = qrshift (Q, R, i, j);
2112%! assert (norm (vec (Q'*Q - eye (3,"single")), Inf), single (0), ...
2113%! 1e1 * eps ("single"));
2114%! assert (norm (vec (triu (R) - R), Inf), single (0));
2115%! assert (norm (vec (Q*R - AA(:,p)), Inf), single (0), ...
2116%! norm (AA)*1e1 * eps ("single"));
2117*/
2118
2119OCTAVE_END_NAMESPACE(octave)
ComplexColumnVector conj(const ComplexColumnVector &a)
Template for N-dimensional array classes with like-type math operators.
Definition MArray.h:61
octave_idx_type length() const
Definition ovl.h:111
bool isinteger() const
Definition ov.h:730
SparseMatrix sparse_matrix_value(bool frc_str_conv=false) const
Definition ov.h:909
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:877
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
bool issparse() const
Definition ov.h:753
FloatMatrix float_matrix_value(bool frc_str_conv=false) const
Definition ov.h:862
bool iscomplex() const
Definition ov.h:741
Array< octave_idx_type > octave_idx_type_vector_value(bool req_int=false, bool frc_str_conv=false, bool frc_vec_conv=false) const
octave_idx_type columns() const
Definition ov.h:547
Matrix matrix_value(bool frc_str_conv=false) const
Definition ov.h:859
FloatComplexMatrix float_complex_matrix_value(bool frc_str_conv=false) const
Definition ov.h:881
SparseComplexMatrix sparse_complex_matrix_value(bool frc_str_conv=false) const
Definition ov.h:913
Definition qr.h:39
OCTAVE_BEGIN_NAMESPACE(octave) static octave_value daspk_fcn
void print_usage()
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:1003
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
octave_value_list ovl(const OV_Args &... args)
Construct an octave_value_list with less typing.
Definition ovl.h:217