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