GNU Octave 11.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-2026 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 "builtin-defun-decls.h"
39#include "defun.h"
40#include "error.h"
41#include "errwarn.h"
42#include "ov.h"
43#include "ovl.h"
44
46
47/*
48## Restore all rand* "state" values
49%!function restore_rand_states (state)
50%! rand ("state", state.rand);
51%! randn ("state", state.randn);
52%!endfunction
53
54%!shared old_state, restore_state
55%! ## Save and restore the states of both random number generators that are
56%! ## tested by the unit tests in this file.
57%! old_state.rand = rand ("state");
58%! old_state.randn = randn ("state");
59%! restore_state = onCleanup (@() restore_rand_states (old_state));
60*/
61
62template <typename MT>
63static octave_value
64get_qr_r (const math::qr<MT>& fact)
65{
66 MT R = fact.R ();
67 if (R.issquare () && fact.regular ())
69 else
70 return R;
71}
72
73template <typename T>
74static typename math::qr<T>::type
75qr_type (int nargout, bool economy)
76{
77 if (nargout == 0 || nargout == 1)
78 return math::qr<T>::raw;
79 else if (economy)
80 return math::qr<T>::economy;
81 else
82 return math::qr<T>::std;
83}
84
85// dense X
86//
87// [Q, R] = qr (X): form Q unitary and R upper triangular such
88// that Q * R = X
89//
90// [Q, R] = qr (X, 0): form the economy decomposition such that if X is
91// m by n then only the first n columns of Q are
92// computed.
93//
94// [Q, R, P] = qr (X): form QRP factorization of X where
95// P is a permutation matrix such that
96// A * P = Q * R
97//
98// [Q, R, P] = qr (X, 0): form the economy decomposition with
99// permutation vector P such that Q * R = X(:, P)
100//
101// qr (X) alone returns the output of the LAPACK routine dgeqrf, such
102// that R = triu (qr (X))
103//
104// sparse X
105//
106// X = qr (A, B): if M < N, X is the minimum 2-norm solution of
107// A\B. If M >= N, X is the least squares
108// approximation of A\B. X is calculated by
109// SPQR-function SuiteSparseQR_min2norm.
110//
111DEFUN (qr, args, nargout,
112 doc: /* -*- texinfo -*-
113@deftypefn {} {[@var{Q}, @var{R}] =} qr (@var{A})
114@deftypefnx {} {[@var{Q}, @var{R}, @var{P}] =} qr (@var{A})
115@deftypefnx {} {@var{R} =} qr (@var{A})
116@deftypefnx {} {[@var{C}, @var{R}] =} qr (@var{A}, @var{B})
117@deftypefnx {} {[@var{C}, @var{R}, @var{P}] =} qr (@var{A}, @var{B})
118@deftypefnx {} {@var{X} =} qr (@var{A}, @var{B}) # sparse only
119@deftypefnx {} {[@dots{}] =} qr (@dots{}, "econ")
120@deftypefnx {} {[@dots{}] =} qr (@dots{}, "vector")
121@deftypefnx {} {[@dots{}] =} qr (@dots{}, "matrix")
122@deftypefnx {} {[@dots{}] =} qr (@dots{}, 0)
123@cindex QR factorization
124Compute the QR@tie{}factorization of @var{A}.
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 @xresult{}
170 1.0000 2.0000
171 3.0000 4.0000
172@end group
173@end example
174
175If just @strong{one return value} is requested then it is @var{R}. (Note:
176unlike most commands, the single return value is not the first return value when
177multiple values are requested.)
178
179If a third output @var{P} is requested, then @code{qr} calculates the permuted
180QR@tie{}factorization
181@tex
182$QR = AP$ where $Q$ is an orthogonal matrix, $R$ is upper triangular, and $P$
183is a permutation matrix.
184@end tex
185@ifnottex
186
187@example
188@var{Q} * @var{R} = @var{A} * @var{P}
189@end example
190
191@noindent
192where @var{Q} is an orthogonal matrix, @var{R} is upper triangular, and
193@var{P} is a permutation matrix.
194@end ifnottex
195
196If @var{A} is dense, the QR@tie{}factorization is computed using standard
197@sc{lapack} subroutines. In this case, the permuted factorization has the
198additional property that the diagonal entries of @var{R} are ordered by
199decreasing magnitude. In other words, @code{abs (diag (@var{R}))} will be
200ordered from largest to smallest.
201
202If @var{A} is sparse, @var{P} is a fill-reducing ordering of the columns of
203@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 @w{@code{@var{R} = chol (@var{A}' * @var{A})}}.
240
241If @var{A} is dense, and an additional input matrix @var{B} is supplied, and
242two return values are requested, then @code{qr} returns @var{C}, where
243@w{@code{@var{C} = @var{Q}' * @var{B}}}. This allows the least squares
244approximation of @w{@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 @w{M < N}, then @var{X} is the minimum 2-norm solution of
256@w{@code{@var{A} \ @var{B}}}. If @w{M >= N}, @var{X} is the least squares
257approximation of @w{@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 matrix.
265In 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 @w{M > N}, then the economy factorization will calculate just
277N rows in @var{R} and N columns in @var{Q} and omit the zeros in @var{R}. If
278@w{M @leq{} N}, there is no difference between the economy and standard
279factorizations.
280
281If the optional argument is the numeric value @code{0} then @code{qr} acts as
282if the @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@w{@code{[@var{Q}, @var{R}, @var{P}] = qr (@var{A})}} allows the construction
311of an orthogonal 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 = Ftriu (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 = Ftriu (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 // Note: LAPACK routine returns matrix with R in
631 // upper-triangular portion. The lower-triangular portion
632 // contains information about Q, but is not relevant here.
633 // Matlab returns just upper-triangular R since v2022a.
634 // Octave will do the same. See bug #64861.
635 retval = Ftriu (ovl (fact.R ()));
636 }
637 break;
638
639 case 2:
640 {
641 math::qr<Matrix> fact (m, type);
642 retval = ovl (fact.Q (), get_qr_r (fact));
643 if (have_b)
644 {
645 if (is_cmplx)
646 retval(0) = fact.Q ().transpose ()
647 * args(1).complex_matrix_value ();
648 else
649 retval(0) = fact.Q ().transpose ()
650 * args(1).matrix_value ();
651 }
652 }
653 break;
654
655 default:
656 {
657 math::qrp<Matrix> fact (m, type);
658 if (vector_p)
659 retval = ovl (fact.Q (), get_qr_r (fact), fact.Pvec ());
660 else
661 retval = ovl (fact.Q (), get_qr_r (fact), fact.P ());
662 }
663 break;
664 }
665 }
666 else if (arg.iscomplex ())
667 {
668 math::qr<ComplexMatrix>::type type
669 = qr_type<ComplexMatrix> (nargout, economy);
670
672
673 switch (nargout)
674 {
675 case 0:
676 case 1:
677 {
678 math::qr<ComplexMatrix> fact (m, type);
679 retval = Ftriu (ovl (fact.R ()));
680 }
681 break;
682
683 case 2:
684 {
685 math::qr<ComplexMatrix> fact (m, type);
686 retval = ovl (fact.Q (), get_qr_r (fact));
687 if (have_b)
688 retval(0) = conj (fact.Q ().transpose ())
689 * args(1).complex_matrix_value ();
690 }
691 break;
692
693 default:
694 {
695 math::qrp<ComplexMatrix> fact (m, type);
696 if (vector_p)
697 retval = ovl (fact.Q (), get_qr_r (fact), fact.Pvec ());
698 else
699 retval = ovl (fact.Q (), get_qr_r (fact), fact.P ());
700 }
701 break;
702 }
703 }
704 else
705 err_wrong_type_arg ("qr", arg);
706 }
707 }
708
709 return retval;
710}
711
712/*
713%!test
714%! a = [0, 2, 1; 2, 1, 2];
715%!
716%! [q, r] = qr (a);
717%! [qe, re] = qr (a, "econ");
718%! [qe2, re2] = qr (a, 0);
719%!
720%! assert (q * r, a, sqrt (eps));
721%! assert (qe * re, a, sqrt (eps));
722%! assert (qe2 * re2, a, sqrt (eps));
723
724%!test
725%! a = [0, 2, 1; 2, 1, 2];
726%!
727%! [q, r, p] = qr (a); # FIXME: not giving right dimensions.
728%! [qe, re, pe] = qr (a, "econ");
729%! [qe2, re2, pe2] = qr (a, 0);
730%!
731%! assert (q * r, a * p, sqrt (eps));
732%! assert (qe * re, a * pe, sqrt (eps));
733%! assert (qe2 * re2, a(:, pe2), sqrt (eps));
734
735%!test
736%! a = [0, 2; 2, 1; 1, 2];
737%!
738%! [q, r] = qr (a);
739%! [qe, re] = qr (a, "econ");
740%! [qe2, re2] = qr (a, 0);
741%!
742%! assert (q * r, a, sqrt (eps));
743%! assert (qe * re, a, sqrt (eps));
744%! assert (qe2 * re2, a, sqrt (eps));
745
746%!test
747%! a = [0, 2; 2, 1; 1, 2];
748%!
749%! [q, r, p] = qr (a);
750%! [qe, re, pe] = qr (a, "econ");
751%! [qe2, re2, pe2] = qr (a, 0);
752%!
753%! assert (q * r, a * p, sqrt (eps));
754%! assert (qe * re, a * pe, sqrt (eps));
755%! assert (qe2 * re2, a(:, pe2), sqrt (eps));
756
757%!test
758%! a = [0, 2, 1; 2, 1, 2; 3, 1, 2];
759%! b = [1, 3, 2; 1, 1, 0; 3, 0, 2];
760%!
761%! [q, r] = qr (a);
762%! [c, re] = qr (a, b);
763%!
764%! assert (r, re, sqrt (eps));
765%! assert (q'*b, c, sqrt (eps));
766
767%!test
768%! a = [0, 2, i; 2, 1, 2; 3, 1, 2];
769%! b = [1, 3, 2; 1, i, 0; 3, 0, 2];
770%!
771%! [q, r] = qr (a);
772%! [c, re] = qr (a, b);
773%!
774%! assert (r, re, sqrt (eps));
775%! assert (q'*b, c, sqrt (eps));
776
777%!test
778%! a = [0, 2, i; 2, 1, 2; 3, 1, 2];
779%! b = [1, 3, 2; 1, 1, 0; 3, 0, 2];
780%!
781%! [q, r] = qr (a);
782%! [c, re] = qr (a, b);
783%!
784%! assert (r, re, sqrt (eps));
785%! assert (q'*b, c, sqrt (eps));
786
787%!test
788%! a = [0, 2, 1; 2, 1, 2; 3, 1, 2];
789%! b = [1, 3, 2; 1, i, 0; 3, 0, 2];
790%!
791%! [q, r] = qr (a);
792%! [c, re] = qr (a, b);
793%!
794%! assert (r, re, sqrt (eps));
795%! assert (q'*b, c, sqrt (eps));
796
797## Empty matrices
798%!test
799%! assert (qr (zeros (0, 0)), zeros (0, 0));
800%! assert (qr (zeros (1, 0)), zeros (1, 0));
801%! assert (qr (zeros (0, 1)), zeros (0, 1));
802
803%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE") <*63069>
804%! assert (qr (sparse (0, 0)), sparse (0, 0));
805%! assert (qr (sparse (1, 0)), sparse (1, 0));
806%! assert (qr (sparse (0, 1)), sparse (0, 1));
807
808%!test <*66488>
809%! ## Orientation of 'p' output for dense matrices
810%! [q, r, p] = qr (eye (3));
811%! assert (size (p), [3, 3]);
812%! [q, r, p] = qr (eye (3), 'vector');
813%! assert (size (p), [1, 3]);
814%! [q, r, p] = qr (eye (3), 'econ');
815%! assert (size (p), [3, 3]);
816%! [q, r, p] = qr (eye (3), 0);
817%! assert (size (p), [1, 3]);
818
819%!test <*64861>
820%! ## Return upper-triangular R for dense matrices
821%! A = [1; 4; 8];
822%! R = qr (A);
823%! assert (R, [-9; 0; 0]);
824%! A = single ([1; 4; 8]);
825%! R = qr (A);
826%! assert (R, single ([-9; 0; 0]));
827%! A = [1i; 4; 8];
828%! R = qr (A);
829%! assert (R, [-9; 0; 0]);
830%! A = single ([1i; 4; 8]);
831%! R = qr (A);
832%! assert (R, single ([-9; 0; 0]));
833
834%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE") <*66488>
835%! ## Orientation of 'p' output for sparse matrices
836%! [q, r, p] = qr (speye (3));
837%! assert (size (p), [3, 3]);
838%! [q, r, p] = qr (speye (3), 'vector');
839%! assert (size (p), [1, 3]);
840
841%!testif HAVE_SPQR, HAVE_CHOLMOD
842%! [q, r, p] = qr (speye (3), 'econ');
843%! assert (size (p), [3, 3]);
844%! [q, r, p] = qr (speye (3), 0);
845%! assert (size (p), [1, 3]);
846
847## Test input validation
848%!error <Invalid call> qr ()
849%!error <Invalid call> qr (1,2,3,4)
850%!error <too many output arguments> [a,b,c,d] = qr (1)
851%!error <option string must be .*, not "foo"> qr (magic (3), "foo")
852%!error <option string must be .*, not "foo"> qr (magic (3), ones (3, 1), "foo")
853%!error <too many output arguments when called with A and B>
854%! [q, r, p] = qr (ones (3, 2), ones (3, 1));
855%!error <too many output arguments when called with A and B>
856%! [q, r, p] = qr (ones (3, 2), ones (3, 1), 0);
857
858%!function retval = __testqr (q, r, a, p)
859%! tol = 100* eps (class (q));
860%! retval = 0;
861%! if (nargin == 3)
862%! n1 = norm (q*r - a);
863%! n2 = norm (q'*q - eye (columns (q)));
864%! retval = (n1 < tol && n2 < tol);
865%! else
866%! n1 = norm (q'*q - eye (columns (q)));
867%! retval = (n1 < tol);
868%! if (isvector (p))
869%! n2 = norm (q*r - a(:,p));
870%! retval = (retval && n2 < tol);
871%! else
872%! n2 = norm (q*r - a*p);
873%! retval = (retval && n2 < tol);
874%! endif
875%! endif
876%!endfunction
877
878%!test
879%! t = ones (24, 1);
880%! j = 1;
881%!
882%! if (false) # eliminate big matrix tests
883%! a = rand (5000, 20);
884%! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
885%! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
886%! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
887%! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
888%!
889%! a = a+1i*eps;
890%! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
891%! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
892%! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
893%! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
894%! endif
895%!
896%! a = [ ones(1,15); sqrt(eps)*eye(15) ];
897%! [q,r] = qr (a); t(j++) = __testqr (q, r, a);
898%! [q,r] = qr (a'); t(j++) = __testqr (q, r, a');
899%! [q,r,p] = qr (a); t(j++) = __testqr (q, r, a, p);
900%! [q,r,p] = qr (a'); t(j++) = __testqr (q, r, a', p);
901%!
902%! a = a+1i*eps;
903%! [q,r] = qr (a); t(j++) = __testqr (q, r, a);
904%! [q,r] = qr (a'); t(j++) = __testqr (q, r, a');
905%! [q,r,p] = qr (a); t(j++) = __testqr (q, r, a, p);
906%! [q,r,p] = qr (a'); t(j++) = __testqr (q, r, a', p);
907%!
908%! a = [ ones(1,15); sqrt(eps)*eye(15) ];
909%! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
910%! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
911%! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
912%! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
913%!
914%! a = a+1i*eps;
915%! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
916%! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
917%! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
918%! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
919%!
920%! a = [ 611 196 -192 407 -8 -52 -49 29
921%! 196 899 113 -192 -71 -43 -8 -44
922%! -192 113 899 196 61 49 8 52
923%! 407 -192 196 611 8 44 59 -23
924%! -8 -71 61 8 411 -599 208 208
925%! -52 -43 49 44 -599 411 208 208
926%! -49 -8 8 59 208 208 99 -911
927%! 29 -44 52 -23 208 208 -911 99 ];
928%! [q,r] = qr (a);
929%!
930%! assert (all (t));
931%! if (all (t))
932%! assert (norm (q*r - a), 0, eps (sum (size (a)) * norm (a)));
933%! endif
934
935%!test
936%! a = single ([0, 2, 1; 2, 1, 2]);
937%!
938%! [q, r] = qr (a);
939%! [qe, re] = qr (a, 0);
940%!
941%! assert (q * r, a, sqrt (eps ("single")));
942%! assert (qe * re, a, sqrt (eps ("single")));
943
944%!test
945%! a = single ([0, 2, 1; 2, 1, 2]);
946%!
947%! [q, r, p] = qr (a); # FIXME: not giving right dimensions.
948%! [qe, re, pe] = qr (a, 0);
949%!
950%! assert (q * r, a * p, sqrt (eps ("single")));
951%! assert (qe * re, a(:, pe), sqrt (eps ("single")));
952
953%!test
954%! a = single ([0, 2; 2, 1; 1, 2]);
955%!
956%! [q, r] = qr (a);
957%! [qe, re] = qr (a, 0);
958%!
959%! assert (q * r, a, sqrt (eps ("single")));
960%! assert (qe * re, a, sqrt (eps ("single")));
961
962%!test
963%! a = single ([0, 2; 2, 1; 1, 2]);
964%!
965%! [q, r, p] = qr (a);
966%! [qe, re, pe] = qr (a, 0);
967%!
968%! assert (q * r, a * p, sqrt (eps ("single")));
969%! assert (qe * re, a(:, pe), sqrt (eps ("single")));
970
971%!test
972%! a = single ([0, 2, 1; 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 ("single")));
979%! assert (q'*b, c, sqrt (eps ("single")));
980
981%!test
982%! a = single ([0, 2, i; 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%! a = single ([0, 2, i; 2, 1, 2; 3, 1, 2]);
993%! b = single ([1, 3, 2; 1, 1, 0; 3, 0, 2]);
994%!
995%! [q, r] = qr (a);
996%! [c, re] = qr (a, b);
997%!
998%! assert (r, re, sqrt (eps));
999%! assert (q'*b, c, sqrt (eps));
1000
1001%!test
1002%! a = single ([0, 2, 1; 2, 1, 2; 3, 1, 2]);
1003%! b = single ([1, 3, 2; 1, i, 0; 3, 0, 2]);
1004%!
1005%! [q, r] = qr (a);
1006%! [c, re] = qr (a, b);
1007%!
1008%! assert (r, re, sqrt (eps ("single")));
1009%! assert (q'*b, c, sqrt (eps ("single")));
1010
1011%!test
1012%! t = ones (24, 1);
1013%! j = 1;
1014%!
1015%! if (false) # eliminate big matrix tests
1016%! a = rand (5000,20);
1017%! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
1018%! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
1019%! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
1020%! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
1021%!
1022%! a = a+1i* eps ("single");
1023%! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
1024%! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
1025%! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
1026%! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
1027%! endif
1028%!
1029%! a = [ ones(1,15); sqrt(eps("single"))*eye(15) ];
1030%! [q,r] = qr (a); t(j++) = __testqr (q, r, a);
1031%! [q,r] = qr (a'); t(j++) = __testqr (q, r, a');
1032%! [q,r,p] = qr (a); t(j++) = __testqr (q, r, a, p);
1033%! [q,r,p] = qr (a'); t(j++) = __testqr (q, r, a', p);
1034%!
1035%! a = a+1i* eps ("single");
1036%! [q,r] = qr (a); t(j++) = __testqr (q, r, a);
1037%! [q,r] = qr (a'); t(j++) = __testqr (q, r, a');
1038%! [q,r,p] = qr (a); t(j++) = __testqr (q, r, a, p);
1039%! [q,r,p] = qr (a'); t(j++) = __testqr (q, r, a', p);
1040%!
1041%! a = [ ones(1,15); sqrt(eps("single"))*eye(15) ];
1042%! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
1043%! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
1044%! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
1045%! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a', p);
1046%!
1047%! a = a+1i* eps ("single");
1048%! [q,r] = qr (a, 0); t(j++) = __testqr (q, r, a);
1049%! [q,r] = qr (a',0); t(j++) = __testqr (q, r, a');
1050%! [q,r,p] = qr (a, 0); t(j++) = __testqr (q, r, a, p);
1051%! [q,r,p] = qr (a',0); t(j++) = __testqr (q, r, a',p);
1052%!
1053%! a = [ 611 196 -192 407 -8 -52 -49 29
1054%! 196 899 113 -192 -71 -43 -8 -44
1055%! -192 113 899 196 61 49 8 52
1056%! 407 -192 196 611 8 44 59 -23
1057%! -8 -71 61 8 411 -599 208 208
1058%! -52 -43 49 44 -599 411 208 208
1059%! -49 -8 8 59 208 208 99 -911
1060%! 29 -44 52 -23 208 208 -911 99 ];
1061%! [q,r] = qr (a);
1062%!
1063%! assert (all (t));
1064%! if (all (t))
1065%! assert (norm (q*r - a), 0, 5000 * eps ("single"));
1066%! endif
1067
1068## The deactivated tests below can't be tested till rectangular back-subs is
1069## implemented for sparse matrices.
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%! r = qr (a);
1078%! assert (r'*r, a'*a, 1e-10);
1079
1080%!testif HAVE_COLAMD; (__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%! q = symamd (a);
1087%! a = a(q,q);
1088%! r = qr (a);
1089%! assert (r'*r, a'*a, 1e-10);
1090
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,d) + speye (n,n);
1097%! [c, r] = qr (a, ones (n,1));
1098%! assert (r\c, full (a)\ones (n,1), 10e-10);
1099
1100%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1101%! n = 20; d = 0.2;
1102%! ## initialize generators to make behavior reproducible
1103%! rand ("state", 42);
1104%! randn ("state", 42);
1105%! a = sprandn (n,n,d) + speye (n,n);
1106%! b = randn (n,2);
1107%! [c, r] = qr (a, b);
1108%! assert (r\c, full (a)\b, 10e-10);
1109
1110## Test under-determined systems!!
1111%!#testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1112%! n = 20; d = 0.2;
1113%! ## initialize generators to make behavior reproducible
1114%! rand ("state", 42);
1115%! randn ("state", 42);
1116%! a = sprandn (n,n+1,d) + speye (n,n+1);
1117%! b = randn (n,2);
1118%! [c, r] = qr (a, b);
1119%! assert (r\c, full (a)\b, 10e-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%! r = qr (a);
1128%! assert (r'*r, a'*a, 1e-10);
1129
1130%!testif HAVE_COLAMD; (__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%! q = symamd (a);
1137%! a = a(q,q);
1138%! r = qr (a);
1139%! assert (r'*r, a'*a, 1e-10);
1140
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,d) + speye (n,n);
1147%! [c, r] = qr (a, ones (n,1));
1148%! assert (r\c, full (a)\ones (n,1), 10e-10);
1149
1150%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1151%! n = 20; d = 0.2;
1152%! ## initialize generators to make behavior reproducible
1153%! rand ("state", 42);
1154%! randn ("state", 42);
1155%! a = 1i* sprandn (n,n,d) + speye (n,n);
1156%! b = randn (n,2);
1157%! [c, r] = qr (a, b);
1158%! assert (r\c, full (a)\b, 10e-10);
1159
1160## Test under-determined systems!!
1161%!#testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1162%! n = 20; d = 0.2;
1163%! ## initialize generators to make behavior reproducible
1164%! rand ("state", 42);
1165%! randn ("state", 42);
1166%! a = 1i* sprandn (n,n+1,d) + speye (n,n+1);
1167%! b = randn (n, 2);
1168%! [c, r] = qr (a, b);
1169%! assert (r\c, full (a)\b, 10e-10);
1170
1171%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1172%! n = 12; m = 20; d = 0.2;
1173%! ## initialize generators to make behavior reproducible
1174%! rand ("state", 42);
1175%! randn ("state", 42);
1176%! a = sprandn (m, n, d);
1177%! b = randn (m, 2);
1178%! [c, r] = qr (a, b);
1179%! assert (r\c, full (a)\b, 10e-10);
1180
1181%!testif HAVE_SPQR, HAVE_CHOLMOD
1182%! n = 12; m = 20; d = 0.2;
1183%! ## initialize generators to make behavior reproducible
1184%! rand ("state", 42);
1185%! randn ("state", 42);
1186%! a = sprandn (m, n, d);
1187%! b = sprandn (m, 2, d);
1188%! [c, r] = qr (a, b, 0);
1189%! [c2, r2] = qr (full (a), full (b), 0);
1190%! assert (r\c, r2\c2, 10e-10);
1191
1192%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1193%! n = 12; m = 20; d = 0.2;
1194%! ## initialize generators to make behavior reproducible
1195%! rand ("state", 42);
1196%! randn ("state", 42);
1197%! a = sprandn (m, n, d);
1198%! b = randn (m, 2);
1199%! [c, r, p] = qr (a, b, "matrix");
1200%! x = p * (r\c);
1201%! [c2, r2] = qr (full (a), b);
1202%! x2 = r2 \ c2;
1203%! assert (x, x2, 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%! [q, r, p] = qr (a, "matrix");
1212%! assert (q * r, a * p, 10e-10);
1213
1214%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
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 = sprandn (m, n, d);
1220%! b = randn (m, 2);
1221%! x = qr (a, b);
1222%! [c2, r2] = qr (full (a), b);
1223%! assert (x, r2\c2, 10e-10);
1224
1225%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
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 = sprandn (m, n, d);
1231%! b = i * randn (m, 2);
1232%! x = qr (a, b);
1233%! [c2, r2] = qr (full (a), b);
1234%! assert (x, r2\c2, 10e-10);
1235
1236%!#testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1237%! n = 12; m = 20; d = 0.2;
1238%! ## initialize generators to make behavior reproducible
1239%! rand ("state", 42);
1240%! randn ("state", 42);
1241%! a = sprandn (m, n, d);
1242%! b = i * randn (m, 2);
1243%! [c, r] = qr (a, b);
1244%! [c2, r2] = qr (full (a), b);
1245%! assert (r\c, r2\c2, 10e-10);
1246
1247%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
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 = sprandn (m, n, d);
1253%! b = i * randn (m, 2);
1254%! [c, r, p] = qr (a, b, "matrix");
1255%! x = p * (r\c);
1256%! [c2, r2] = qr (full (a), b);
1257%! x2 = r2 \ c2;
1258%! assert (x, x2, 10e-10);
1259
1260%!testif HAVE_SPQR, HAVE_CHOLMOD
1261%! n = 12; m = 20; d = 0.2;
1262%! ## initialize generators to make behavior reproducible
1263%! rand ("state", 42);
1264%! randn ("state", 42);
1265%! a = i * sprandn (m, n, d);
1266%! b = sprandn (m, 2, d);
1267%! [c, r] = qr (a, b, 0);
1268%! [c2, r2] = qr (full (a), full (b), 0);
1269%! assert (r\c, r2\c2, 10e-10);
1270
1271%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1272%! n = 12; m = 20; d = 0.2;
1273%! ## initialize generators to make behavior reproducible
1274%! rand ("state", 42);
1275%! randn ("state", 42);
1276%! a = i * sprandn (m, n, d);
1277%! b = randn (m, 2);
1278%! [c, r, p] = qr (a, b, "matrix");
1279%! x = p * (r\c);
1280%! [c2, r2] = qr (full (a), b);
1281%! x2 = r2 \ c2;
1282%! assert(x, x2, 10e-10);
1283
1284%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1285%! n = 12; m = 20; d = 0.2;
1286%! ## initialize generators to make behavior reproducible
1287%! rand ("state", 42);
1288%! randn ("state", 42);
1289%! a = i * sprandn (m, n, d);
1290%! [q, r, p] = qr (a, "matrix");
1291%! assert(q * r, a * p, 10e-10);
1292
1293%!testif ; (__have_feature__ ("SPQR") && __have_feature__ ("CHOLMOD")) || __have_feature__ ("CXSPARSE")
1294%! n = 12; m = 20; d = 0.2;
1295%! ## initialize generators to make behavior reproducible
1296%! rand ("state", 42);
1297%! randn ("state", 42);
1298%! a = i * sprandn (m, n, d);
1299%! b = randn (m, 2);
1300%! x = qr (a, b);
1301%! [c2, r2] = qr (full (a), b);
1302%! assert (x, r2\c2, 10e-10);
1303
1304%!testif HAVE_SPQR, HAVE_CHOLMOD
1305%! a = sparse (5, 6);
1306%! a(3,1) = 0.8;
1307%! a(2,2) = 1.4;
1308%! a(1,6) = -0.5;
1309%! r = qr (a);
1310%! assert (r'*r, a'*a, 10e-10);
1311*/
1312
1313static
1314bool
1315check_qr_dims (const octave_value& q, const octave_value& r,
1316 bool allow_ecf = false)
1317{
1318 octave_idx_type m = q.rows ();
1319 octave_idx_type k = r.rows ();
1320 octave_idx_type n = r.columns ();
1321 return ((q.ndims () == 2 && r.ndims () == 2 && k == q.columns ())
1322 && (m == k || (allow_ecf && k == n && k < m)));
1323}
1324
1325static
1326bool
1327check_index (const octave_value& i, bool vector_allowed = false)
1328{
1329 return ((i.isreal () || i.isinteger ())
1330 && (i.is_scalar_type () || vector_allowed));
1331}
1332
1333DEFUN (qrupdate, args, ,
1334 doc: /* -*- texinfo -*-
1335@deftypefn {} {[@var{Q1}, @var{R1}] =} qrupdate (@var{Q}, @var{R}, @var{u}, @var{v})
1336Update a QR factorization given update vectors or matrices.
1337
1338Given a QR@tie{}factorization of a real or complex matrix
1339@w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and
1340@var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization of
1341@w{@var{A} + @var{u}*@var{v}'}, where @var{u} and @var{v} are column vectors
1342(rank-1 update) or matrices with equal number of columns
1343(rank-k update). Notice that the latter case is done as a sequence of
1344rank-1 updates; thus, for k large enough, it will be both faster and more
1345accurate to recompute the factorization from scratch.
1346
1347The QR@tie{}factorization supplied may be either full (Q is square) or
1348economized (R is square).
1349
1350@seealso{qr, qrinsert, qrdelete, qrshift}
1351@end deftypefn */)
1352{
1353 octave_value_list retval;
1354
1355 if (args.length () != 4)
1356 print_usage ();
1357
1358 octave_value argq = args(0);
1359 octave_value argr = args(1);
1360 octave_value argu = args(2);
1361 octave_value argv = args(3);
1362
1363 if (! argq.isnumeric () || ! argr.isnumeric ()
1364 || ! argu.isnumeric () || ! argv.isnumeric ())
1365 print_usage ();
1366
1367 if (! check_qr_dims (argq, argr, true))
1368 error ("qrupdate: Q and R dimensions don't match");
1369
1370 if (argq.isreal () && argr.isreal () && argu.isreal ()
1371 && argv.isreal ())
1372 {
1373 // all real case
1374 if (argq.is_single_type () || argr.is_single_type ()
1375 || argu.is_single_type () || argv.is_single_type ())
1376 {
1378 FloatMatrix R = argr.float_matrix_value ();
1379 FloatMatrix u = argu.float_matrix_value ();
1380 FloatMatrix v = argv.float_matrix_value ();
1381
1382 math::qr<FloatMatrix> fact (Q, R);
1383 fact.update (u, v);
1384
1385 retval = ovl (fact.Q (), get_qr_r (fact));
1386 }
1387 else
1388 {
1389 Matrix Q = argq.matrix_value ();
1390 Matrix R = argr.matrix_value ();
1391 Matrix u = argu.matrix_value ();
1392 Matrix v = argv.matrix_value ();
1393
1394 math::qr<Matrix> fact (Q, R);
1395 fact.update (u, v);
1396
1397 retval = ovl (fact.Q (), get_qr_r (fact));
1398 }
1399 }
1400 else
1401 {
1402 // complex case
1403 if (argq.is_single_type () || argr.is_single_type ()
1404 || argu.is_single_type () || argv.is_single_type ())
1405 {
1410
1411 math::qr<FloatComplexMatrix> fact (Q, R);
1412 fact.update (u, v);
1413
1414 retval = ovl (fact.Q (), get_qr_r (fact));
1415 }
1416 else
1417 {
1422
1423 math::qr<ComplexMatrix> fact (Q, R);
1424 fact.update (u, v);
1425
1426 retval = ovl (fact.Q (), get_qr_r (fact));
1427 }
1428 }
1429
1430 return retval;
1431}
1432
1433/*
1434%!shared A, u, v, Ac, uc, vc
1435%! A = [0.091364 0.613038 0.999083;
1436%! 0.594638 0.425302 0.603537;
1437%! 0.383594 0.291238 0.085574;
1438%! 0.265712 0.268003 0.238409;
1439%! 0.669966 0.743851 0.445057 ];
1440%!
1441%! u = [0.85082;
1442%! 0.76426;
1443%! 0.42883;
1444%! 0.53010;
1445%! 0.80683 ];
1446%!
1447%! v = [0.98810;
1448%! 0.24295;
1449%! 0.43167 ];
1450%!
1451%! Ac = [0.620405 + 0.956953i 0.480013 + 0.048806i 0.402627 + 0.338171i;
1452%! 0.589077 + 0.658457i 0.013205 + 0.279323i 0.229284 + 0.721929i;
1453%! 0.092758 + 0.345687i 0.928679 + 0.241052i 0.764536 + 0.832406i;
1454%! 0.912098 + 0.721024i 0.049018 + 0.269452i 0.730029 + 0.796517i;
1455%! 0.112849 + 0.603871i 0.486352 + 0.142337i 0.355646 + 0.151496i ];
1456%!
1457%! uc = [0.20351 + 0.05401i;
1458%! 0.13141 + 0.43708i;
1459%! 0.29808 + 0.08789i;
1460%! 0.69821 + 0.38844i;
1461%! 0.74871 + 0.25821i ];
1462%!
1463%! vc = [0.85839 + 0.29468i;
1464%! 0.20820 + 0.93090i;
1465%! 0.86184 + 0.34689i ];
1466%!
1467
1468%!test
1469%! [Q,R] = qr (A);
1470%! [Q,R] = qrupdate (Q, R, u, v);
1471%! assert (norm (vec (Q'*Q - eye (5)), Inf), 0, 1e1*eps);
1472%! assert (norm (vec (triu (R)-R), Inf), 0);
1473%! assert (norm (vec (Q*R - A - u*v'), Inf), 0, norm (A)*1e1*eps);
1474%!
1475%!test
1476%! [Q,R] = qr (Ac);
1477%! [Q,R] = qrupdate (Q, R, uc, vc);
1478%! assert (norm (vec (Q'*Q - eye (5)), Inf), 0, 1e1*eps);
1479%! assert (norm (vec (triu (R)-R), Inf), 0);
1480%! assert (norm (vec (Q*R - Ac - uc*vc'), Inf), 0, norm (Ac)*1e1*eps);
1481
1482%!test
1483%! [Q,R] = qr (single (A));
1484%! [Q,R] = qrupdate (Q, R, single (u), single (v));
1485%! assert (norm (vec (Q'*Q - eye (5,"single")), Inf), single (0), ...
1486%! 1e1 * eps ("single"));
1487%! assert (norm (vec (triu (R)-R), Inf), single (0));
1488%! assert (norm (vec (Q*R - single (A) - single (u)* single (v)'), Inf), ...
1489%! single (0), norm (single (A))*1e1 * eps ("single"));
1490%!
1491%!test
1492%! [Q,R] = qr (single (Ac));
1493%! [Q,R] = qrupdate (Q, R, single (uc), single (vc));
1494%! assert (norm (vec (Q'*Q - eye (5,"single")), Inf), single (0), ...
1495%! 1e1 * eps ("single"));
1496%! assert (norm (vec (triu (R)-R), Inf), single (0));
1497%! assert (norm (vec (Q*R - single (Ac) - single (uc)* single (vc)'), Inf), ...
1498%! single (0), norm (single (Ac))*1e1 * eps ("single"));
1499*/
1500
1501DEFUN (qrinsert, args, ,
1502 doc: /* -*- texinfo -*-
1503@deftypefn {} {[@var{Q1}, @var{R1}] =} qrinsert (@var{Q}, @var{R}, @var{j}, @var{x}, @var{orient})
1504Update a QR factorization given a row or column to insert in the original
1505factored matrix.
1506
1507
1508Given a QR@tie{}factorization of a real or complex matrix
1509@w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and
1510@var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization of
1511@w{[A(:,1:j-1) x A(:,j:n)]}, where @var{u} is a column vector to be inserted
1512into @var{A} (if @var{orient} is @qcode{"col"}), or the
1513QR@tie{}factorization of @w{[A(1:j-1,:);x;A(:,j:n)]}, where @var{x} is a row
1514vector to be inserted into @var{A} (if @var{orient} is @qcode{"row"}).
1515
1516The default value of @var{orient} is @qcode{"col"}. If @var{orient} is
1517@qcode{"col"}, @var{u} may be a matrix and @var{j} an index vector
1518resulting in the QR@tie{}factorization of a matrix @var{B} such that
1519@w{B(:,@var{j})}@ gives @var{u} and @w{B(:,@var{j}) = []}@ gives @var{A}.
1520Notice that the latter case is done as a sequence of k insertions;
1521thus, for k large enough, it will be both faster and more accurate to
1522recompute the factorization from scratch.
1523
1524If @var{orient} is @qcode{"col"}, the QR@tie{}factorization supplied may
1525be either full (Q is square) or economized (R is square).
1526
1527If @var{orient} is @qcode{"row"}, full factorization is needed.
1528@seealso{qr, qrupdate, qrdelete, qrshift}
1529@end deftypefn */)
1530{
1531 int nargin = args.length ();
1532
1533 if (nargin < 4 || nargin > 5)
1534 print_usage ();
1535
1536 octave_value argq = args(0);
1537 octave_value argr = args(1);
1538 octave_value argj = args(2);
1539 octave_value argx = args(3);
1540
1541 if (! argq.isnumeric () || ! argr.isnumeric ()
1542 || ! argx.isnumeric ()
1543 || (nargin > 4 && ! args(4).is_string ()))
1544 print_usage ();
1545
1546 std::string orient = (nargin < 5) ? "col" : args(4).string_value ();
1547 bool col = (orient == "col");
1548
1549 if (! col && orient != "row")
1550 error (R"(qrinsert: ORIENT must be "col" or "row")");
1551
1552 if (! check_qr_dims (argq, argr, col) || (! col && argx.rows () != 1))
1553 error ("qrinsert: dimension mismatch");
1554
1555 if (! check_index (argj, col))
1556 error ("qrinsert: invalid index J");
1557
1558 octave_value_list retval;
1559
1561
1562 octave_idx_type one = 1;
1563
1564 if (argq.isreal () && argr.isreal () && argx.isreal ())
1565 {
1566 // real case
1567 if (argq.is_single_type () || argr.is_single_type ()
1568 || argx.is_single_type ())
1569 {
1571 FloatMatrix R = argr.float_matrix_value ();
1573
1574 math::qr<FloatMatrix> fact (Q, R);
1575
1576 if (col)
1577 fact.insert_col (x, j-one);
1578 else
1579 fact.insert_row (x.row (0), j(0)-one);
1580
1581 retval = ovl (fact.Q (), get_qr_r (fact));
1582 }
1583 else
1584 {
1585 Matrix Q = argq.matrix_value ();
1586 Matrix R = argr.matrix_value ();
1587 Matrix x = argx.matrix_value ();
1588
1589 math::qr<Matrix> 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 }
1599 else
1600 {
1601 // complex case
1602 if (argq.is_single_type () || argr.is_single_type ()
1603 || argx.is_single_type ())
1604 {
1608
1609 math::qr<FloatComplexMatrix> fact (Q, R);
1610
1611 if (col)
1612 fact.insert_col (x, j-one);
1613 else
1614 fact.insert_row (x.row (0), j(0)-one);
1615
1616 retval = ovl (fact.Q (), get_qr_r (fact));
1617 }
1618 else
1619 {
1623
1624 math::qr<ComplexMatrix> fact (Q, R);
1625
1626 if (col)
1627 fact.insert_col (x, j-one);
1628 else
1629 fact.insert_row (x.row (0), j(0)-one);
1630
1631 retval = ovl (fact.Q (), get_qr_r (fact));
1632 }
1633 }
1634
1635 return retval;
1636}
1637
1638/*
1639%!test
1640%! [Q,R] = qr (A);
1641%! [Q,R] = qrinsert (Q, R, 3, u);
1642%! assert (norm (vec (Q'*Q - eye (5)), Inf), 0, 1e1*eps);
1643%! assert (norm (vec (triu (R) - R), Inf), 0);
1644%! assert (norm (vec (Q*R - [A(:,1:2) u A(:,3)]), Inf), 0, norm (A)*1e1*eps);
1645%!test
1646%! [Q,R] = qr (Ac);
1647%! [Q,R] = qrinsert (Q, R, 3, uc);
1648%! assert (norm (vec (Q'*Q - eye (5)), Inf), 0, 1e1*eps);
1649%! assert (norm (vec (triu (R) - R), Inf), 0);
1650%! assert (norm (vec (Q*R - [Ac(:,1:2) uc Ac(:,3)]), Inf), 0, ...
1651%! norm (Ac) * 1e1 * eps);
1652%!test
1653%! x = [0.85082 0.76426 0.42883 ];
1654%!
1655%! [Q,R] = qr (A);
1656%! [Q,R] = qrinsert (Q, R, 3, x, "row");
1657%! assert (norm (vec (Q'*Q - eye (6)), Inf), 0, 1e1*eps);
1658%! assert (norm (vec (triu (R) - R), Inf), 0);
1659%! assert (norm (vec (Q*R - [A(1:2,:);x;A(3:5,:)]), Inf), 0, norm (A)*1e1*eps);
1660%!test
1661%! x = [0.20351 + 0.05401i 0.13141 + 0.43708i 0.29808 + 0.08789i ];
1662%!
1663%! [Q,R] = qr (Ac);
1664%! [Q,R] = qrinsert (Q, R, 3, x, "row");
1665%! assert (norm (vec (Q'*Q - eye (6)), Inf), 0, 1e1*eps);
1666%! assert (norm (vec (triu (R) - R), Inf), 0);
1667%! assert (norm (vec (Q*R - [Ac(1:2,:);x;Ac(3:5,:)]), Inf), 0, ...
1668%! norm (Ac) * 1e1 * eps);
1669
1670%!test
1671%! [Q,R] = qr (single (A));
1672%! [Q,R] = qrinsert (Q, R, 3, single (u));
1673%! assert (norm (vec (Q'*Q - eye (5,"single")), Inf), single (0), ...
1674%! 1e1 * eps ("single"));
1675%! assert (norm (vec (triu (R) - R), Inf), single (0));
1676%! assert (norm (vec (Q*R - single ([A(:,1:2) u A(:,3)])), Inf), ...
1677%! single (0), norm (single (A))*1e1 * eps ("single"));
1678%!test
1679%! [Q,R] = qr (single (Ac));
1680%! [Q,R] = qrinsert (Q, R, 3, single (uc));
1681%! assert (norm (vec (Q'*Q - eye (5,"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) uc Ac(:,3)])), Inf), ...
1685%! single (0), norm (single (Ac))*1e1 * eps ("single"));
1686%!test
1687%! x = single ([0.85082 0.76426 0.42883 ]);
1688%!
1689%! [Q,R] = qr (single (A));
1690%! [Q,R] = qrinsert (Q, R, 3, x, "row");
1691%! assert (norm (vec (Q'*Q - eye (6,"single")), Inf), single (0), ...
1692%! 1e1 * eps ("single"));
1693%! assert (norm (vec (triu (R) - R), Inf), single (0));
1694%! assert (norm (vec (Q*R - single ([A(1:2,:);x;A(3:5,:)])), Inf), ...
1695%! single (0), norm (single (A))*1e1 * eps ("single"));
1696%!test
1697%! x = single ([0.20351 + 0.05401i 0.13141 + 0.43708i 0.29808 + 0.08789i ]);
1698%!
1699%! [Q,R] = qr (single (Ac));
1700%! [Q,R] = qrinsert (Q, R, 3, x, "row");
1701%! assert (norm (vec (Q'*Q - eye (6,"single")), Inf), single (0), ...
1702%! 1e1 * eps ("single"));
1703%! assert (norm (vec (triu (R) - R), Inf), single (0));
1704%! assert (norm (vec (Q*R - single ([Ac(1:2,:);x;Ac(3:5,:)])), Inf), ...
1705%! single (0), norm (single (Ac))*1e1 * eps ("single"));
1706*/
1707
1708DEFUN (qrdelete, args, ,
1709 doc: /* -*- texinfo -*-
1710@deftypefn {} {[@var{Q1}, @var{R1}] =} qrdelete (@var{Q}, @var{R}, @var{j}, @var{orient})
1711Update a QR factorization given a row or column to delete from the original
1712factored matrix.
1713
1714Given a QR@tie{}factorization of a real or complex matrix
1715@w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and
1716@var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization of
1717@w{[A(:,1:j-1), U, A(:,j:n)]},
1718where @var{u} is a column vector to be inserted into @var{A}
1719(if @var{orient} is @qcode{"col"}),
1720or the QR@tie{}factorization of @w{[A(1:j-1,:);X;A(:,j:n)]},
1721where @var{x} is a row @var{orient} is @qcode{"row"}).
1722The default value of @var{orient} is @qcode{"col"}.
1723
1724If @var{orient} is @qcode{"col"}, @var{j} may be an index vector
1725resulting in the QR@tie{}factorization of a matrix @var{B} such that
1726@w{A(:,@var{j}) = []}@ gives @var{B}. Notice that the latter case is done as
1727a sequence of k deletions; thus, for k large enough, it will be both faster
1728and more accurate to recompute the factorization from scratch.
1729
1730If @var{orient} is @qcode{"col"}, the QR@tie{}factorization supplied may
1731be either full (Q is square) or economized (R is square).
1732
1733If @var{orient} is @qcode{"row"}, full factorization is needed.
1734@seealso{qr, qrupdate, qrinsert, qrshift}
1735@end deftypefn */)
1736{
1737 int nargin = args.length ();
1738
1739 if (nargin < 3 || nargin > 4)
1740 print_usage ();
1741
1742 octave_value argq = args(0);
1743 octave_value argr = args(1);
1744 octave_value argj = args(2);
1745
1746 if (! argq.isnumeric () || ! argr.isnumeric ()
1747 || (nargin > 3 && ! args(3).is_string ()))
1748 print_usage ();
1749
1750 std::string orient = (nargin < 4) ? "col" : args(3).string_value ();
1751 bool col = orient == "col";
1752
1753 if (! col && orient != "row")
1754 error (R"(qrdelete: ORIENT must be "col" or "row")");
1755
1756 if (! check_qr_dims (argq, argr, col))
1757 error ("qrdelete: dimension mismatch");
1758
1760 if (! check_index (argj, col))
1761 error ("qrdelete: invalid index J");
1762
1763 octave_value_list retval;
1764
1765 octave_idx_type one = 1;
1766
1767 if (argq.isreal () && argr.isreal ())
1768 {
1769 // real case
1770 if (argq.is_single_type () || argr.is_single_type ())
1771 {
1773 FloatMatrix R = argr.float_matrix_value ();
1774
1775 math::qr<FloatMatrix> 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 else
1785 {
1786 Matrix Q = argq.matrix_value ();
1787 Matrix R = argr.matrix_value ();
1788
1789 math::qr<Matrix> fact (Q, R);
1790
1791 if (col)
1792 fact.delete_col (j-one);
1793 else
1794 fact.delete_row (j(0)-one);
1795
1796 retval = ovl (fact.Q (), get_qr_r (fact));
1797 }
1798 }
1799 else
1800 {
1801 // complex case
1802 if (argq.is_single_type () || argr.is_single_type ())
1803 {
1806
1807 math::qr<FloatComplexMatrix> fact (Q, R);
1808
1809 if (col)
1810 fact.delete_col (j-one);
1811 else
1812 fact.delete_row (j(0)-one);
1813
1814 retval = ovl (fact.Q (), get_qr_r (fact));
1815 }
1816 else
1817 {
1820
1821 math::qr<ComplexMatrix> fact (Q, R);
1822
1823 if (col)
1824 fact.delete_col (j-one);
1825 else
1826 fact.delete_row (j(0)-one);
1827
1828 retval = ovl (fact.Q (), get_qr_r (fact));
1829 }
1830 }
1831
1832 return retval;
1833}
1834
1835/*
1836%!test
1837%! AA = [0.091364 0.613038 0.027504 0.999083;
1838%! 0.594638 0.425302 0.562834 0.603537;
1839%! 0.383594 0.291238 0.742073 0.085574;
1840%! 0.265712 0.268003 0.783553 0.238409;
1841%! 0.669966 0.743851 0.457255 0.445057 ];
1842%!
1843%! [Q,R] = qr (AA);
1844%! [Q,R] = qrdelete (Q, R, 3);
1845%! assert (norm (vec (Q'*Q - eye (5)), Inf), 0, 16*eps);
1846%! assert (norm (vec (triu (R) - R), Inf), 0);
1847%! assert (norm (vec (Q*R - [AA(:,1:2) AA(:,4)]), Inf), 0, norm (AA)*1e1*eps);
1848%!
1849%!test
1850%! AA = [0.364554 + 0.993117i 0.669818 + 0.510234i 0.426568 + 0.041337i 0.847051 + 0.233291i;
1851%! 0.049600 + 0.242783i 0.448946 + 0.484022i 0.141155 + 0.074420i 0.446746 + 0.392706i;
1852%! 0.581922 + 0.657416i 0.581460 + 0.030016i 0.219909 + 0.447288i 0.201144 + 0.069132i;
1853%! 0.694986 + 0.000571i 0.682327 + 0.841712i 0.807537 + 0.166086i 0.192767 + 0.358098i;
1854%! 0.945002 + 0.066788i 0.350492 + 0.642638i 0.579629 + 0.048102i 0.600170 + 0.636938i ] * I;
1855%!
1856%! [Q,R] = qr (AA);
1857%! [Q,R] = qrdelete (Q, R, 3);
1858%! assert (norm (vec (Q'*Q - eye (5)), Inf) < 16*eps);
1859%! assert (norm (vec (triu (R) - R), Inf) == 0);
1860%! assert (norm (vec (Q*R - [AA(:,1:2) AA(:,4)]), Inf) < norm (AA)*1e1*eps);
1861%!
1862%!test
1863%! AA = [0.091364 0.613038 0.027504 0.999083;
1864%! 0.594638 0.425302 0.562834 0.603537;
1865%! 0.383594 0.291238 0.742073 0.085574;
1866%! 0.265712 0.268003 0.783553 0.238409;
1867%! 0.669966 0.743851 0.457255 0.445057 ];
1868%!
1869%! [Q,R] = qr (AA);
1870%! [Q,R] = qrdelete (Q, R, 3, "row");
1871%! assert (norm (vec (Q'*Q - eye (4)), Inf), 0, 1e1*eps);
1872%! assert (norm (vec (triu (R) - R), Inf), 0);
1873%! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf), 0, norm (AA)*1e1*eps);
1874%!
1875%!test
1876%! AA = [0.364554 + 0.993117i 0.669818 + 0.510234i 0.426568 + 0.041337i 0.847051 + 0.233291i;
1877%! 0.049600 + 0.242783i 0.448946 + 0.484022i 0.141155 + 0.074420i 0.446746 + 0.392706i;
1878%! 0.581922 + 0.657416i 0.581460 + 0.030016i 0.219909 + 0.447288i 0.201144 + 0.069132i;
1879%! 0.694986 + 0.000571i 0.682327 + 0.841712i 0.807537 + 0.166086i 0.192767 + 0.358098i;
1880%! 0.945002 + 0.066788i 0.350492 + 0.642638i 0.579629 + 0.048102i 0.600170 + 0.636938i ] * I;
1881%!
1882%! [Q,R] = qr (AA);
1883%! [Q,R] = qrdelete (Q, R, 3, "row");
1884%! assert (norm (vec (Q'*Q - eye (4)), Inf), 0, 1e1*eps);
1885%! assert (norm (vec (triu (R) - R), Inf), 0);
1886%! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf), 0, norm (AA)*1e1*eps);
1887
1888%!test
1889%! AA = single ([0.091364 0.613038 0.027504 0.999083;
1890%! 0.594638 0.425302 0.562834 0.603537;
1891%! 0.383594 0.291238 0.742073 0.085574;
1892%! 0.265712 0.268003 0.783553 0.238409;
1893%! 0.669966 0.743851 0.457255 0.445057 ]);
1894%!
1895%! [Q,R] = qr (AA);
1896%! [Q,R] = qrdelete (Q, R, 3);
1897%! assert (norm (vec (Q'*Q - eye (5,"single")), Inf), single (0), ...
1898%! 1e1 * eps ("single"));
1899%! assert (norm (vec (triu (R) - R), Inf), single (0));
1900%! assert (norm (vec (Q*R - [AA(:,1:2) AA(:,4)]), Inf), single (0), ...
1901%! norm (AA)*1e1 * eps ("single"));
1902%!
1903%!test
1904%! AA = single ([0.364554 + 0.993117i 0.669818 + 0.510234i 0.426568 + 0.041337i 0.847051 + 0.233291i;
1905%! 0.049600 + 0.242783i 0.448946 + 0.484022i 0.141155 + 0.074420i 0.446746 + 0.392706i;
1906%! 0.581922 + 0.657416i 0.581460 + 0.030016i 0.219909 + 0.447288i 0.201144 + 0.069132i;
1907%! 0.694986 + 0.000571i 0.682327 + 0.841712i 0.807537 + 0.166086i 0.192767 + 0.358098i;
1908%! 0.945002 + 0.066788i 0.350492 + 0.642638i 0.579629 + 0.048102i 0.600170 + 0.636938i ]) * I;
1909%!
1910%! [Q,R] = qr (AA);
1911%! [Q,R] = qrdelete (Q, R, 3);
1912%! assert (norm (vec (Q'*Q - eye (5,"single")), Inf), single (0), ...
1913%! 1e1 * eps ("single"));
1914%! assert (norm (vec (triu (R) - R), Inf), single (0));
1915%! assert (norm (vec (Q*R - [AA(:,1:2) AA(:,4)]), Inf), single (0), ...
1916%! norm (AA)*1e1 * eps ("single"));
1917
1918%!test
1919%! AA = single ([0.091364 0.613038 0.027504 0.999083;
1920%! 0.594638 0.425302 0.562834 0.603537;
1921%! 0.383594 0.291238 0.742073 0.085574;
1922%! 0.265712 0.268003 0.783553 0.238409;
1923%! 0.669966 0.743851 0.457255 0.445057 ]);
1924%!
1925%! [Q,R] = qr (AA);
1926%! [Q,R] = qrdelete (Q, R, 3, "row");
1927%! assert (norm (vec (Q'*Q - eye (4,"single")), Inf), single (0), ...
1928%! 1.5e1 * eps ("single"));
1929%! assert (norm (vec (triu (R) - R), Inf), single (0));
1930%! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf), single (0), ...
1931%! norm (AA)*1e1 * eps ("single"));
1932%!testif HAVE_QRUPDATE
1933%! ## Same test as above but with more precicision
1934%! AA = single ([0.091364 0.613038 0.027504 0.999083;
1935%! 0.594638 0.425302 0.562834 0.603537;
1936%! 0.383594 0.291238 0.742073 0.085574;
1937%! 0.265712 0.268003 0.783553 0.238409;
1938%! 0.669966 0.743851 0.457255 0.445057 ]);
1939%!
1940%! [Q,R] = qr (AA);
1941%! [Q,R] = qrdelete (Q, R, 3, "row");
1942%! assert (norm (vec (Q'*Q - eye (4,"single")), Inf), single (0), ...
1943%! 1e1* eps ("single"));
1944%! assert (norm (vec (triu (R) - R), Inf), single (0));
1945%! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf), single (0), ...
1946%! norm (AA)*1e1 * eps ("single"));
1947%!
1948%!test
1949%! AA = single ([0.364554 + 0.993117i 0.669818 + 0.510234i 0.426568 + 0.041337i 0.847051 + 0.233291i;
1950%! 0.049600 + 0.242783i 0.448946 + 0.484022i 0.141155 + 0.074420i 0.446746 + 0.392706i;
1951%! 0.581922 + 0.657416i 0.581460 + 0.030016i 0.219909 + 0.447288i 0.201144 + 0.069132i;
1952%! 0.694986 + 0.000571i 0.682327 + 0.841712i 0.807537 + 0.166086i 0.192767 + 0.358098i;
1953%! 0.945002 + 0.066788i 0.350492 + 0.642638i 0.579629 + 0.048102i 0.600170 + 0.636938i ]) * I;
1954%!
1955%! [Q,R] = qr (AA);
1956%! [Q,R] = qrdelete (Q, R, 3, "row");
1957%! assert (norm (vec (Q'*Q - eye (4,"single")), Inf), single (0), ...
1958%! 1e1 * eps ("single"));
1959%! assert (norm (vec (triu (R) - R), Inf), single (0));
1960%! assert (norm (vec (Q*R - [AA(1:2,:);AA(4:5,:)]), Inf), single (0), ...
1961%! norm (AA)*1e1 * eps ("single"));
1962*/
1963
1964DEFUN (qrshift, args, ,
1965 doc: /* -*- texinfo -*-
1966@deftypefn {} {[@var{Q1}, @var{R1}] =} qrshift (@var{Q}, @var{R}, @var{i}, @var{j})
1967Update a QR factorization given a range of columns to shift in the original
1968factored matrix.
1969
1970Given a QR@tie{}factorization of a real or complex matrix
1971@w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and
1972@var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization
1973of @w{@var{A}(:,p)}, where @w{p} is the permutation @*
1974@code{p = [1:i-1, circshift(i:j, 1), j+1:n]} if @w{@var{i} < @var{j}} @*
1975 or @*
1976@code{p = [1:j-1, circshift(j:i,-1), i+1:n]} if @w{@var{j} < @var{i}}. @*
1977
1978@seealso{qr, qrupdate, qrinsert, qrdelete}
1979@end deftypefn */)
1980{
1981 if (args.length () != 4)
1982 print_usage ();
1983
1984 octave_value argq = args(0);
1985 octave_value argr = args(1);
1986 octave_value argi = args(2);
1987 octave_value argj = args(3);
1988
1989 if (! argq.isnumeric () || ! argr.isnumeric ())
1990 print_usage ();
1991
1992 if (! check_qr_dims (argq, argr, true))
1993 error ("qrshift: dimensions mismatch");
1994
1995 octave_idx_type i = argi.idx_type_value ();
1996 octave_idx_type j = argj.idx_type_value ();
1997
1998 if (! check_index (argi) || ! check_index (argj))
1999 error ("qrshift: invalid index I or J");
2000
2001 octave_value_list retval;
2002
2003 if (argq.isreal () && argr.isreal ())
2004 {
2005 // all real case
2006 if (argq.is_single_type ()
2007 && argr.is_single_type ())
2008 {
2010 FloatMatrix R = argr.float_matrix_value ();
2011
2012 math::qr<FloatMatrix> fact (Q, R);
2013 fact.shift_cols (i-1, j-1);
2014
2015 retval = ovl (fact.Q (), get_qr_r (fact));
2016 }
2017 else
2018 {
2019 Matrix Q = argq.matrix_value ();
2020 Matrix R = argr.matrix_value ();
2021
2022 math::qr<Matrix> fact (Q, R);
2023 fact.shift_cols (i-1, j-1);
2024
2025 retval = ovl (fact.Q (), get_qr_r (fact));
2026 }
2027 }
2028 else
2029 {
2030 // complex case
2031 if (argq.is_single_type ()
2032 && argr.is_single_type ())
2033 {
2036
2037 math::qr<FloatComplexMatrix> fact (Q, R);
2038 fact.shift_cols (i-1, j-1);
2039
2040 retval = ovl (fact.Q (), get_qr_r (fact));
2041 }
2042 else
2043 {
2046
2047 math::qr<ComplexMatrix> fact (Q, R);
2048 fact.shift_cols (i-1, j-1);
2049
2050 retval = ovl (fact.Q (), get_qr_r (fact));
2051 }
2052 }
2053
2054 return retval;
2055}
2056
2057/*
2058%!test
2059%! AA = A.';
2060%! i = 2; j = 4; p = [1:i-1, circshift(i:j,-1), j+1:5];
2061%!
2062%! [Q,R] = qr (AA);
2063%! [Q,R] = qrshift (Q, R, i, j);
2064%! assert (norm (vec (Q'*Q - eye (3)), Inf), 0, 1e1*eps);
2065%! assert (norm (vec (triu (R) - R), Inf), 0);
2066%! assert (norm (vec (Q*R - AA(:,p)), Inf), 0, norm (AA)*1e1*eps);
2067%!
2068%! j = 2; i = 4; p = [1:j-1, circshift(j:i,+1), i+1:5];
2069%!
2070%! [Q,R] = qr (AA);
2071%! [Q,R] = qrshift (Q, R, i, j);
2072%! assert (norm (vec (Q'*Q - eye (3)), Inf), 0, 1e1*eps);
2073%! assert (norm (vec (triu (R) - R), Inf), 0);
2074%! assert (norm (vec (Q*R - AA(:,p)), Inf), 0, norm (AA)*1e1*eps);
2075%!
2076%!test
2077%! AA = Ac.';
2078%! i = 2; j = 4; p = [1:i-1, circshift(i:j,-1), j+1:5];
2079%!
2080%! [Q,R] = qr (AA);
2081%! [Q,R] = qrshift (Q, R, i, j);
2082%! assert (norm (vec (Q'*Q - eye (3)), Inf), 0, 1e1*eps);
2083%! assert (norm (vec (triu (R) - R), Inf), 0);
2084%! assert (norm (vec (Q*R - AA(:,p)), Inf), 0, norm (AA)*1e1*eps);
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)), Inf), 0, 1e1*eps);
2091%! assert (norm (vec (triu (R) - R), Inf), 0);
2092%! assert (norm (vec (Q*R - AA(:,p)), Inf), 0, norm (AA)*1e1*eps);
2093
2094%!test
2095%! AA = single (A).';
2096%! i = 2; j = 4; p = [1:i-1, circshift(i:j,-1), j+1:5];
2097%!
2098%! [Q,R] = qr (AA);
2099%! [Q,R] = qrshift (Q, R, i, j);
2100%! assert (norm (vec (Q'*Q - eye (3,"single")), Inf), single (0), ...
2101%! 1e1 * eps ("single"));
2102%! assert (norm (vec (triu (R) - R), Inf), single (0));
2103%! assert (norm (vec (Q*R - AA(:,p)), Inf), single (0), ...
2104%! norm (AA)*1e1 * eps ("single"));
2105%!
2106%! j = 2; i = 4; p = [1:j-1, circshift(j:i,+1), i+1:5];
2107%!
2108%! [Q,R] = qr (AA);
2109%! [Q,R] = qrshift (Q, R, i, j);
2110%! assert (norm (vec (Q'*Q - eye (3,"single")), Inf), single (0), ...
2111%! 1e1 * eps ("single"));
2112%! assert (norm (vec (triu (R) - R), Inf), single (0));
2113%! assert (norm (vec (Q*R - AA(:,p)), Inf), single (0), ...
2114%! norm (AA)*1e1 * eps ("single"));
2115%!
2116%!test
2117%! AA = single (Ac).';
2118%! i = 2; j = 4; p = [1:i-1, circshift(i:j,-1), j+1:5];
2119%!
2120%! [Q,R] = qr (AA);
2121%! [Q,R] = qrshift (Q, R, i, j);
2122%! assert (norm (vec (Q'*Q - eye (3,"single")), Inf), single (0), ...
2123%! 1e1 * eps ("single"));
2124%! assert (norm (vec (triu (R) - R), Inf), single (0));
2125%! assert (norm (vec (Q*R - AA(:,p)), Inf), single (0), ...
2126%! norm (AA)*1e1 * eps ("single"));
2127%!
2128%! j = 2; i = 4; p = [1:j-1, circshift(j:i,+1), i+1:5];
2129%!
2130%! [Q,R] = qr (AA);
2131%! [Q,R] = qrshift (Q, R, i, j);
2132%! assert (norm (vec (Q'*Q - eye (3,"single")), Inf), single (0), ...
2133%! 1e1 * eps ("single"));
2134%! assert (norm (vec (triu (R) - R), Inf), single (0));
2135%! assert (norm (vec (Q*R - AA(:,p)), Inf), single (0), ...
2136%! norm (AA)*1e1 * eps ("single"));
2137*/
2138
2139OCTAVE_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:728
SparseMatrix sparse_matrix_value(bool frc_str_conv=false) const
Definition ov.h:907
int ndims() const
Definition ov.h:549
octave_idx_type rows() const
Definition ov.h:543
bool is_scalar_type() const
Definition ov.h:742
bool isreal() const
Definition ov.h:736
bool isnumeric() const
Definition ov.h:748
ComplexMatrix complex_matrix_value(bool frc_str_conv=false) const
Definition ov.h:875
octave_idx_type idx_type_value(bool req_int=false, bool frc_str_conv=false) const
bool is_single_type() const
Definition ov.h:696
bool issparse() const
Definition ov.h:751
FloatMatrix float_matrix_value(bool frc_str_conv=false) const
Definition ov.h:860
bool iscomplex() const
Definition ov.h:739
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:545
Matrix matrix_value(bool frc_str_conv=false) const
Definition ov.h:857
FloatComplexMatrix float_complex_matrix_value(bool frc_str_conv=false) const
Definition ov.h:879
SparseComplexMatrix sparse_complex_matrix_value(bool frc_str_conv=false) const
Definition ov.h:911
Definition qr.h:38
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:1008
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
octave_value_list ovl(const OV_Args &... args)
Construct an octave_value_list with less typing.
Definition ovl.h:217
F77_RET_T const F77_DBLE * x
octave_value_list Ftriu(const octave_value_list &args, int)
Definition tril.cc:431