qr.cc

Go to the documentation of this file.
00001 /*
00002 
00003 Copyright (C) 1996-2012 John W. Eaton
00004 Copyright (C) 2008-2009 Jaroslav Hajek
00005 Copyright (C) 2008-2009 VZLU Prague
00006 
00007 This file is part of Octave.
00008 
00009 Octave is free software; you can redistribute it and/or modify it
00010 under the terms of the GNU General Public License as published by the
00011 Free Software Foundation; either version 3 of the License, or (at your
00012 option) any later version.
00013 
00014 Octave is distributed in the hope that it will be useful, but WITHOUT
00015 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00016 FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00017 for more details.
00018 
00019 You should have received a copy of the GNU General Public License
00020 along with Octave; see the file COPYING.  If not, see
00021 <http://www.gnu.org/licenses/>.
00022 
00023 */
00024 
00025 #ifdef HAVE_CONFIG_H
00026 #include <config.h>
00027 #endif
00028 
00029 #include "CmplxQR.h"
00030 #include "CmplxQRP.h"
00031 #include "dbleQR.h"
00032 #include "dbleQRP.h"
00033 #include "fCmplxQR.h"
00034 #include "fCmplxQRP.h"
00035 #include "floatQR.h"
00036 #include "floatQRP.h"
00037 #include "SparseQR.h"
00038 #include "SparseCmplxQR.h"
00039 
00040 
00041 #include "defun-dld.h"
00042 #include "error.h"
00043 #include "gripes.h"
00044 #include "oct-obj.h"
00045 #include "utils.h"
00046 
00047 template <class MT>
00048 static octave_value
00049 get_qr_r (const base_qr<MT>& fact)
00050 {
00051   MT R = fact.R ();
00052   if (R.is_square () && fact.regular ())
00053     return octave_value (R, MatrixType (MatrixType::Upper));
00054   else
00055     return R;
00056 }
00057 
00058 // [Q, R] = qr (X):      form Q unitary and R upper triangular such
00059 //                        that Q * R = X
00060 //
00061 // [Q, R] = qr (X, 0):    form the economy decomposition such that if X is
00062 //                        m by n then only the first n columns of Q are
00063 //                        computed.
00064 //
00065 // [Q, R, P] = qr (X):    form QRP factorization of X where
00066 //                        P is a permutation matrix such that
00067 //                        A * P = Q * R
00068 //
00069 // [Q, R, P] = qr (X, 0): form the economy decomposition with
00070 //                        permutation vector P such that Q * R = X (:, P)
00071 //
00072 // qr (X) alone returns the output of the LAPACK routine dgeqrf, such
00073 // that R = triu (qr (X))
00074 
00075 DEFUN_DLD (qr, args, nargout,
00076   "-*- texinfo -*-\n\
00077 @deftypefn  {Loadable Function} {[@var{Q}, @var{R}, @var{P}] =} qr (@var{A})\n\
00078 @deftypefnx {Loadable Function} {[@var{Q}, @var{R}, @var{P}] =} qr (@var{A}, '0')\n\
00079 @deftypefnx {Loadable Function} {[@var{C}, @var{R}] =} qr (@var{A}, @var{B})\n\
00080 @deftypefnx {Loadable Function} {[@var{C}, @var{R}] =} qr (@var{A}, @var{B}, '0')\n\
00081 @cindex QR factorization\n\
00082 Compute the QR@tie{}factorization of @var{A}, using standard @sc{lapack}\n\
00083 subroutines.  For example, given the matrix @code{@var{A} = [1, 2; 3, 4]},\n\
00084 \n\
00085 @example\n\
00086 [@var{Q}, @var{R}] = qr (@var{A})\n\
00087 @end example\n\
00088 \n\
00089 @noindent\n\
00090 returns\n\
00091 \n\
00092 @example\n\
00093 @group\n\
00094 @var{Q} =\n\
00095 \n\
00096   -0.31623  -0.94868\n\
00097   -0.94868   0.31623\n\
00098 \n\
00099 @var{R} =\n\
00100 \n\
00101   -3.16228  -4.42719\n\
00102    0.00000  -0.63246\n\
00103 @end group\n\
00104 @end example\n\
00105 \n\
00106 The @code{qr} factorization has applications in the solution of least\n\
00107 squares problems\n\
00108 @tex\n\
00109 $$\n\
00110 \\min_x \\left\\Vert A x - b \\right\\Vert_2\n\
00111 $$\n\
00112 @end tex\n\
00113 @ifnottex\n\
00114 \n\
00115 @example\n\
00116 @code{min norm(A x - b)}\n\
00117 @end example\n\
00118 \n\
00119 @end ifnottex\n\
00120 for overdetermined systems of equations (i.e.,\n\
00121 @tex\n\
00122 $A$\n\
00123 @end tex\n\
00124 @ifnottex\n\
00125 @var{A}\n\
00126 @end ifnottex\n\
00127  is a tall, thin matrix).  The QR@tie{}factorization is\n\
00128 @tex\n\
00129 $QR = A$ where $Q$ is an orthogonal matrix and $R$ is upper triangular.\n\
00130 @end tex\n\
00131 @ifnottex\n\
00132 @code{@var{Q} * @var{Q} = @var{A}} where @var{Q} is an orthogonal matrix and\n\
00133 @var{R} is upper triangular.\n\
00134 @end ifnottex\n\
00135 \n\
00136 If given a second argument of '0', @code{qr} returns an economy-sized\n\
00137 QR@tie{}factorization, omitting zero rows of @var{R} and the corresponding\n\
00138 columns of @var{Q}.\n\
00139 \n\
00140 If the matrix @var{A} is full, the permuted QR@tie{}factorization\n\
00141 @code{[@var{Q}, @var{R}, @var{P}] = qr (@var{A})} forms the\n\
00142 QR@tie{}factorization such that the diagonal entries of @var{R} are\n\
00143 decreasing in magnitude order.  For example, given the matrix @code{a = [1,\n\
00144 2; 3, 4]},\n\
00145 \n\
00146 @example\n\
00147 [@var{Q}, @var{R}, @var{P}] = qr (@var{A})\n\
00148 @end example\n\
00149 \n\
00150 @noindent\n\
00151 returns\n\
00152 \n\
00153 @example\n\
00154 @group\n\
00155 @var{Q} =\n\
00156 \n\
00157   -0.44721  -0.89443\n\
00158   -0.89443   0.44721\n\
00159 \n\
00160 @var{R} =\n\
00161 \n\
00162   -4.47214  -3.13050\n\
00163    0.00000   0.44721\n\
00164 \n\
00165 @var{P} =\n\
00166 \n\
00167    0  1\n\
00168    1  0\n\
00169 @end group\n\
00170 @end example\n\
00171 \n\
00172 The permuted @code{qr} factorization @code{[@var{Q}, @var{R}, @var{P}] = qr\n\
00173 (@var{A})} factorization allows the construction of an orthogonal basis of\n\
00174 @code{span (A)}.\n\
00175 \n\
00176 If the matrix @var{A} is sparse, then compute the sparse\n\
00177 QR@tie{}factorization of @var{A}, using @sc{CSparse}.  As the matrix @var{Q}\n\
00178 is in general a full matrix, this function returns the @var{Q}-less\n\
00179 factorization @var{R} of @var{A}, such that @code{@var{R} = chol (@var{A}' *\n\
00180 @var{A})}.\n\
00181 \n\
00182 If the final argument is the scalar @code{0} and the number of rows is\n\
00183 larger than the number of columns, then an economy factorization is\n\
00184 returned.  That is @var{R} will have only @code{size (@var{A},1)} rows.\n\
00185 \n\
00186 If an additional matrix @var{B} is supplied, then @code{qr} returns\n\
00187 @var{C}, where @code{@var{C} = @var{Q}' * @var{B}}.  This allows the\n\
00188 least squares approximation of @code{@var{A} \\ @var{B}} to be calculated\n\
00189 as\n\
00190 \n\
00191 @example\n\
00192 @group\n\
00193 [@var{C}, @var{R}] = qr (@var{A}, @var{B})\n\
00194 x = @var{R} \\ @var{C}\n\
00195 @end group\n\
00196 @end example\n\
00197 @end deftypefn")
00198 {
00199   octave_value_list retval;
00200 
00201   int nargin = args.length ();
00202 
00203   if (nargin < 1 || nargin > (args(0).is_sparse_type() ? 3 : 2))
00204     {
00205       print_usage ();
00206       return retval;
00207     }
00208 
00209   octave_value arg = args(0);
00210 
00211   int arg_is_empty = empty_arg ("qr", arg.rows (), arg.columns ());
00212 
00213   if (arg_is_empty < 0)
00214     return retval;
00215 
00216   if (arg.is_sparse_type ())
00217     {
00218       bool economy = false;
00219       bool is_cmplx = false;
00220       int have_b = 0;
00221 
00222       if (arg.is_complex_type ())
00223         is_cmplx = true;
00224       if (nargin > 1)
00225         {
00226           have_b = 1;
00227           if (args(nargin-1).is_scalar_type ())
00228             {
00229               int val = args(nargin-1).int_value ();
00230               if (val == 0)
00231                 {
00232                   economy = true;
00233                   have_b = (nargin > 2 ? 2 : 0);
00234                 }
00235             }
00236           if (have_b > 0 && args(have_b).is_complex_type ())
00237             is_cmplx = true;
00238         }
00239 
00240       if (!error_state)
00241         {
00242           if (have_b && nargout < 2)
00243             error ("qr: incorrect number of output arguments");
00244           else if (is_cmplx)
00245             {
00246               SparseComplexQR q (arg.sparse_complex_matrix_value ());
00247               if (!error_state)
00248                 {
00249                   if (have_b > 0)
00250                     {
00251                       retval(1) = q.R (economy);
00252                       retval(0) = q.C (args(have_b).complex_matrix_value ());
00253                       if (arg.rows() < arg.columns())
00254                         warning ("qr: non minimum norm solution for under-determined problem");
00255                     }
00256                   else if (nargout > 1)
00257                     {
00258                       retval(1) = q.R (economy);
00259                       retval(0) = q.Q ();
00260                     }
00261                   else
00262                     retval(0) = q.R (economy);
00263                 }
00264             }
00265           else
00266             {
00267               SparseQR q (arg.sparse_matrix_value ());
00268               if (!error_state)
00269                 {
00270                   if (have_b > 0)
00271                     {
00272                       retval(1) = q.R (economy);
00273                       retval(0) = q.C (args(have_b).matrix_value ());
00274                       if (args(0).rows() < args(0).columns())
00275                         warning ("qr: non minimum norm solution for under-determined problem");
00276                     }
00277                   else if (nargout > 1)
00278                     {
00279                       retval(1) = q.R (economy);
00280                       retval(0) = q.Q ();
00281                     }
00282                   else
00283                     retval(0) = q.R (economy);
00284                 }
00285             }
00286         }
00287     }
00288   else
00289     {
00290       QR::type type = (nargout == 0 || nargout == 1) ? QR::raw
00291         : (nargin == 2 ? QR::economy : QR::std);
00292 
00293       if (arg.is_single_type ())
00294         {
00295           if (arg.is_real_type ())
00296             {
00297               FloatMatrix m = arg.float_matrix_value ();
00298 
00299               if (! error_state)
00300                 {
00301                   switch (nargout)
00302                     {
00303                     case 0:
00304                     case 1:
00305                       {
00306                         FloatQR fact (m, type);
00307                         retval(0) = fact.R ();
00308                       }
00309                       break;
00310 
00311                     case 2:
00312                       {
00313                         FloatQR fact (m, type);
00314                         retval(1) = get_qr_r (fact);
00315                         retval(0) = fact.Q ();
00316                       }
00317                       break;
00318 
00319                     default:
00320                       {
00321                         FloatQRP fact (m, type);
00322                         if (type == QR::economy)
00323                           retval(2) = fact.Pvec ();
00324                         else
00325                           retval(2) = fact.P ();
00326                         retval(1) = get_qr_r (fact);
00327                         retval(0) = fact.Q ();
00328                       }
00329                       break;
00330                     }
00331                 }
00332             }
00333           else if (arg.is_complex_type ())
00334             {
00335               FloatComplexMatrix m = arg.float_complex_matrix_value ();
00336 
00337               if (! error_state)
00338                 {
00339                   switch (nargout)
00340                     {
00341                     case 0:
00342                     case 1:
00343                       {
00344                         FloatComplexQR fact (m, type);
00345                         retval(0) = fact.R ();
00346                       }
00347                       break;
00348 
00349                     case 2:
00350                       {
00351                         FloatComplexQR fact (m, type);
00352                         retval(1) = get_qr_r (fact);
00353                         retval(0) = fact.Q ();
00354                       }
00355                       break;
00356 
00357                     default:
00358                       {
00359                         FloatComplexQRP fact (m, type);
00360                         if (type == QR::economy)
00361                           retval(2) = fact.Pvec ();
00362                         else
00363                           retval(2) = fact.P ();
00364                         retval(1) = get_qr_r (fact);
00365                         retval(0) = fact.Q ();
00366                       }
00367                       break;
00368                     }
00369                 }
00370             }
00371         }
00372       else
00373         {
00374           if (arg.is_real_type ())
00375             {
00376               Matrix m = arg.matrix_value ();
00377 
00378               if (! error_state)
00379                 {
00380                   switch (nargout)
00381                     {
00382                     case 0:
00383                     case 1:
00384                       {
00385                         QR fact (m, type);
00386                         retval(0) = fact.R ();
00387                       }
00388                       break;
00389 
00390                     case 2:
00391                       {
00392                         QR fact (m, type);
00393                         retval(1) = get_qr_r (fact);
00394                         retval(0) = fact.Q ();
00395                       }
00396                       break;
00397 
00398                     default:
00399                       {
00400                         QRP fact (m, type);
00401                         if (type == QR::economy)
00402                           retval(2) = fact.Pvec ();
00403                         else
00404                           retval(2) = fact.P ();
00405                         retval(1) = get_qr_r (fact);
00406                         retval(0) = fact.Q ();
00407                       }
00408                       break;
00409                     }
00410                 }
00411             }
00412           else if (arg.is_complex_type ())
00413             {
00414               ComplexMatrix m = arg.complex_matrix_value ();
00415 
00416               if (! error_state)
00417                 {
00418                   switch (nargout)
00419                     {
00420                     case 0:
00421                     case 1:
00422                       {
00423                         ComplexQR fact (m, type);
00424                         retval(0) = fact.R ();
00425                       }
00426                       break;
00427 
00428                     case 2:
00429                       {
00430                         ComplexQR fact (m, type);
00431                         retval(1) = get_qr_r (fact);
00432                         retval(0) = fact.Q ();
00433                       }
00434                       break;
00435 
00436                     default:
00437                       {
00438                         ComplexQRP fact (m, type);
00439                         if (type == QR::economy)
00440                           retval(2) = fact.Pvec ();
00441                         else
00442                           retval(2) = fact.P ();
00443                         retval(1) = get_qr_r (fact);
00444                         retval(0) = fact.Q ();
00445                       }
00446                       break;
00447                     }
00448                 }
00449             }
00450           else
00451             gripe_wrong_type_arg ("qr", arg);
00452         }
00453     }
00454 
00455   return retval;
00456 }
00457 
00458 /*
00459 
00460 %!test
00461 %! a = [0, 2, 1; 2, 1, 2];
00462 %!
00463 %! [q, r] = qr (a);
00464 %!
00465 %! [qe, re] = qr (a, 0);
00466 %!
00467 %! assert (q * r, a, sqrt (eps));
00468 %! assert (qe * re, a, sqrt (eps));
00469 
00470 %!test
00471 %! a = [0, 2, 1; 2, 1, 2];
00472 %!
00473 %! [q, r, p] = qr (a);  # not giving right dimensions. FIXME
00474 %!
00475 %! [qe, re, pe] = qr (a, 0);
00476 %!
00477 %! assert (q * r, a * p, sqrt (eps));
00478 %! assert (qe * re, a(:, pe), sqrt (eps));
00479 
00480 %!test
00481 %! a = [0, 2; 2, 1; 1, 2];
00482 %!
00483 %! [q, r] = qr (a);
00484 %!
00485 %! [qe, re] = qr (a, 0);
00486 %!
00487 %! assert (q * r, a, sqrt (eps));
00488 %! assert (qe * re, a, sqrt (eps));
00489 
00490 %!test
00491 %! a = [0, 2; 2, 1; 1, 2];
00492 %!
00493 %! [q, r, p] = qr (a);
00494 %!
00495 %! [qe, re, pe] = qr (a, 0);
00496 %!
00497 %! assert (q * r, a * p, sqrt (eps));
00498 %! assert (qe * re, a(:, pe), sqrt (eps));
00499 
00500 %!error <Invalid call to qr> qr ();
00501 %!error <Invalid call to qr> qr ([1, 2; 3, 4], 0, 2);
00502 
00503 %!function retval = __testqr (q, r, a, p)
00504 %!  tol = 100*eps (class(q));
00505 %!  retval = 0;
00506 %!  if (nargin == 3)
00507 %!    n1 = norm (q*r-a);
00508 %!    n2 = norm (q'*q-eye(columns(q)));
00509 %!    retval = (n1 < tol && n2 < tol);
00510 %!  else
00511 %!    n1 = norm (q'*q-eye(columns(q)));
00512 %!    retval = (n1 < tol);
00513 %!    if (isvector (p))
00514 %!      n2 = norm (q*r-a(:,p));
00515 %!      retval = (retval && n2 < tol);
00516 %!    else
00517 %!      n2 = norm (q*r - a*p);
00518 %!      retval = (retval && n2 < tol);
00519 %!    endif
00520 %!  endif
00521 %!endfunction
00522 
00523 %!test
00524 %!
00525 %! t = ones (24, 1);
00526 %! j = 1;
00527 %!
00528 %! if false # eliminate big matrix tests
00529 %!   a = rand(5000,20);
00530 %!   [q,r]=qr(a,0); t(j++) = __testqr(q,r,a);
00531 %!   [q,r]=qr(a',0); t(j++) = __testqr(q,r,a');
00532 %!   [q,r,p]=qr(a,0); t(j++) = __testqr(q,r,a,p);
00533 %!   [q,r,p]=qr(a',0); t(j++) = __testqr(q,r,a',p);
00534 %!
00535 %!   a = a+1i*eps;
00536 %!   [q,r]=qr(a,0); t(j++) = __testqr(q,r,a);
00537 %!   [q,r]=qr(a',0); t(j++) = __testqr(q,r,a');
00538 %!   [q,r,p]=qr(a,0); t(j++) = __testqr(q,r,a,p);
00539 %!   [q,r,p]=qr(a',0); t(j++) = __testqr(q,r,a',p);
00540 %! endif
00541 %!
00542 %! a = [ ones(1,15); sqrt(eps)*eye(15) ];
00543 %! [q,r]=qr(a); t(j++) = __testqr(q,r,a);
00544 %! [q,r]=qr(a'); t(j++) = __testqr(q,r,a');
00545 %! [q,r,p]=qr(a); t(j++) = __testqr(q,r,a,p);
00546 %! [q,r,p]=qr(a'); t(j++) = __testqr(q,r,a',p);
00547 %!
00548 %! a = a+1i*eps;
00549 %! [q,r]=qr(a); t(j++) = __testqr(q,r,a);
00550 %! [q,r]=qr(a'); t(j++) = __testqr(q,r,a');
00551 %! [q,r,p]=qr(a); t(j++) = __testqr(q,r,a,p);
00552 %! [q,r,p]=qr(a'); t(j++) = __testqr(q,r,a',p);
00553 %!
00554 %! a = [ ones(1,15); sqrt(eps)*eye(15) ];
00555 %! [q,r]=qr(a,0); t(j++) = __testqr(q,r,a);
00556 %! [q,r]=qr(a',0); t(j++) = __testqr(q,r,a');
00557 %! [q,r,p]=qr(a,0); t(j++) = __testqr(q,r,a,p);
00558 %! [q,r,p]=qr(a',0); t(j++) = __testqr(q,r,a',p);
00559 %!
00560 %! a = a+1i*eps;
00561 %! [q,r]=qr(a,0); t(j++) = __testqr(q,r,a);
00562 %! [q,r]=qr(a',0); t(j++) = __testqr(q,r,a');
00563 %! [q,r,p]=qr(a,0); t(j++) = __testqr(q,r,a,p);
00564 %! [q,r,p]=qr(a',0); t(j++) = __testqr(q,r,a',p);
00565 %!
00566 %! a = [
00567 %! 611   196  -192   407    -8   -52   -49    29
00568 %! 196   899   113  -192   -71   -43    -8   -44
00569 %! -192   113   899   196    61    49     8    52
00570 %! 407  -192   196   611     8    44    59   -23
00571 %! -8   -71    61     8   411  -599   208   208
00572 %! -52   -43    49    44  -599   411   208   208
00573 %! -49    -8     8    59   208   208    99  -911
00574 %! 29   -44    52   -23   208   208  -911    99
00575 %! ];
00576 %! [q,r] = qr(a);
00577 %!
00578 %! assert(all (t) && norm(q*r-a) < 5000*eps);
00579 
00580 %!test
00581 %! a = single ([0, 2, 1; 2, 1, 2]);
00582 %!
00583 %! [q, r] = qr (a);
00584 %!
00585 %! [qe, re] = qr (a, 0);
00586 %!
00587 %! assert (q * r, a, sqrt (eps ('single')));
00588 %! assert (qe * re, a, sqrt (eps ('single')));
00589 
00590 %!test
00591 %! a = single([0, 2, 1; 2, 1, 2]);
00592 %!
00593 %! [q, r, p] = qr (a);  # not giving right dimensions. FIXME
00594 %!
00595 %! [qe, re, pe] = qr (a, 0);
00596 %!
00597 %! assert (q * r, a * p, sqrt (eps('single')));
00598 %! assert (qe * re, a(:, pe), sqrt (eps('single')));
00599 
00600 %!test
00601 %! a = single([0, 2; 2, 1; 1, 2]);
00602 %!
00603 %! [q, r] = qr (a);
00604 %!
00605 %! [qe, re] = qr (a, 0);
00606 %!
00607 %! assert (q * r, a, sqrt (eps('single')));
00608 %! assert (qe * re, a, sqrt (eps('single')));
00609 
00610 %!test
00611 %! a = single([0, 2; 2, 1; 1, 2]);
00612 %!
00613 %! [q, r, p] = qr (a);
00614 %!
00615 %! [qe, re, pe] = qr (a, 0);
00616 %!
00617 %! assert (q * r, a * p, sqrt (eps('single')));
00618 %! assert (qe * re, a(:, pe), sqrt (eps('single')));
00619 
00620 %!error <Invalid call to qr> qr ();
00621 %!error <Invalid call to qr> qr ([1, 2; 3, 4], 0, 2);
00622 
00623 %!test
00624 %!
00625 %! t = ones (24, 1);
00626 %! j = 1;
00627 %!
00628 %! if false # eliminate big matrix tests
00629 %!   a = rand(5000,20);
00630 %!   [q,r]=qr(a,0); t(j++) = __testqr(q,r,a);
00631 %!   [q,r]=qr(a',0); t(j++) = __testqr(q,r,a');
00632 %!   [q,r,p]=qr(a,0); t(j++) = __testqr(q,r,a,p);
00633 %!   [q,r,p]=qr(a',0); t(j++) = __testqr(q,r,a',p);
00634 %!
00635 %!   a = a+1i*eps('single');
00636 %!   [q,r]=qr(a,0); t(j++) = __testqr(q,r,a);
00637 %!   [q,r]=qr(a',0); t(j++) = __testqr(q,r,a');
00638 %!   [q,r,p]=qr(a,0); t(j++) = __testqr(q,r,a,p);
00639 %!   [q,r,p]=qr(a',0); t(j++) = __testqr(q,r,a',p);
00640 %! endif
00641 %!
00642 %! a = [ ones(1,15); sqrt(eps('single'))*eye(15) ];
00643 %! [q,r]=qr(a); t(j++) = __testqr(q,r,a);
00644 %! [q,r]=qr(a'); t(j++) = __testqr(q,r,a');
00645 %! [q,r,p]=qr(a); t(j++) = __testqr(q,r,a,p);
00646 %! [q,r,p]=qr(a'); t(j++) = __testqr(q,r,a',p);
00647 %!
00648 %! a = a+1i*eps('single');
00649 %! [q,r]=qr(a); t(j++) = __testqr(q,r,a);
00650 %! [q,r]=qr(a'); t(j++) = __testqr(q,r,a');
00651 %! [q,r,p]=qr(a); t(j++) = __testqr(q,r,a,p);
00652 %! [q,r,p]=qr(a'); t(j++) = __testqr(q,r,a',p);
00653 %!
00654 %! a = [ ones(1,15); sqrt(eps('single'))*eye(15) ];
00655 %! [q,r]=qr(a,0); t(j++) = __testqr(q,r,a);
00656 %! [q,r]=qr(a',0); t(j++) = __testqr(q,r,a');
00657 %! [q,r,p]=qr(a,0); t(j++) = __testqr(q,r,a,p);
00658 %! [q,r,p]=qr(a',0); t(j++) = __testqr(q,r,a',p);
00659 %!
00660 %! a = a+1i*eps('single');
00661 %! [q,r]=qr(a,0); t(j++) = __testqr(q,r,a);
00662 %! [q,r]=qr(a',0); t(j++) = __testqr(q,r,a');
00663 %! [q,r,p]=qr(a,0); t(j++) = __testqr(q,r,a,p);
00664 %! [q,r,p]=qr(a',0); t(j++) = __testqr(q,r,a',p);
00665 %!
00666 %! a = [
00667 %! 611   196  -192   407    -8   -52   -49    29
00668 %! 196   899   113  -192   -71   -43    -8   -44
00669 %! -192   113   899   196    61    49     8    52
00670 %! 407  -192   196   611     8    44    59   -23
00671 %! -8   -71    61     8   411  -599   208   208
00672 %! -52   -43    49    44  -599   411   208   208
00673 %! -49    -8     8    59   208   208    99  -911
00674 %! 29   -44    52   -23   208   208  -911    99
00675 %! ];
00676 %! [q,r] = qr(a);
00677 %!
00678 %! assert(all (t) && norm(q*r-a) < 5000*eps('single'));
00679 
00680 %% The deactivated tests below can't be tested till rectangular back-subs is
00681 %% implemented for sparse matrices.
00682 
00683 %!testif HAVE_CXSPARSE
00684 %! n = 20; d= 0.2;
00685 %! a = sprandn(n,n,d)+speye(n,n);
00686 %! r = qr(a);
00687 %! assert(r'*r,a'*a,1e-10)
00688 
00689 %!testif HAVE_COLAMD
00690 %! n = 20; d= 0.2;
00691 %! a = sprandn(n,n,d)+speye(n,n);
00692 %! q = symamd(a);
00693 %! a = a(q,q);
00694 %! r = qr(a);
00695 %! assert(r'*r,a'*a,1e-10)
00696 
00697 %!testif HAVE_CXSPARSE
00698 %! n = 20; d= 0.2;
00699 %! a = sprandn(n,n,d)+speye(n,n);
00700 %! [c,r] = qr(a,ones(n,1));
00701 %! assert (r\c,full(a)\ones(n,1),10e-10)
00702 
00703 %!testif HAVE_CXSPARSE
00704 %! n = 20; d= 0.2;
00705 %! a = sprandn(n,n,d)+speye(n,n);
00706 %! b = randn(n,2);
00707 %! [c,r] = qr(a,b);
00708 %! assert (r\c,full(a)\b,10e-10)
00709 
00710 %% Test under-determined systems!!
00711 %!#testif HAVE_CXSPARSE
00712 %! n = 20; d= 0.2;
00713 %! a = sprandn(n,n+1,d)+speye(n,n+1);
00714 %! b = randn(n,2);
00715 %! [c,r] = qr(a,b);
00716 %! assert (r\c,full(a)\b,10e-10)
00717 
00718 %!testif HAVE_CXSPARSE
00719 %! n = 20; d= 0.2;
00720 %! a = 1i*sprandn(n,n,d)+speye(n,n);
00721 %! r = qr(a);
00722 %! assert(r'*r,a'*a,1e-10)
00723 
00724 %!testif HAVE_COLAMD
00725 %! n = 20; d= 0.2;
00726 %! a = 1i*sprandn(n,n,d)+speye(n,n);
00727 %! q = symamd(a);
00728 %! a = a(q,q);
00729 %! r = qr(a);
00730 %! assert(r'*r,a'*a,1e-10)
00731 
00732 %!testif HAVE_CXSPARSE
00733 %! n = 20; d= 0.2;
00734 %! a = 1i*sprandn(n,n,d)+speye(n,n);
00735 %! [c,r] = qr(a,ones(n,1));
00736 %! assert (r\c,full(a)\ones(n,1),10e-10)
00737 
00738 %!testif HAVE_CXSPARSE
00739 %! n = 20; d= 0.2;
00740 %! a = 1i*sprandn(n,n,d)+speye(n,n);
00741 %! b = randn(n,2);
00742 %! [c,r] = qr(a,b);
00743 %! assert (r\c,full(a)\b,10e-10)
00744 
00745 %% Test under-determined systems!!
00746 %!#testif HAVE_CXSPARSE
00747 %! n = 20; d= 0.2;
00748 %! a = 1i*sprandn(n,n+1,d)+speye(n,n+1);
00749 %! b = randn(n,2);
00750 %! [c,r] = qr(a,b);
00751 %! assert (r\c,full(a)\b,10e-10)
00752 
00753 %!error qr(sprandn(10,10,0.2),ones(10,1));
00754 
00755 */
00756 
00757 static
00758 bool check_qr_dims (const octave_value& q, const octave_value& r,
00759                     bool allow_ecf = false)
00760 {
00761   octave_idx_type m = q.rows (), k = r.rows (), n = r.columns ();
00762   return ((q.ndims () == 2 && r.ndims () == 2 && k == q.columns ())
00763             && (m == k || (allow_ecf && k == n && k < m)));
00764 }
00765 
00766 static
00767 bool check_index (const octave_value& i, bool vector_allowed = false)
00768 {
00769   return ((i.is_real_type () || i.is_integer_type ())
00770           && (i.is_scalar_type () || vector_allowed));
00771 }
00772 
00773 DEFUN_DLD (qrupdate, args, ,
00774   "-*- texinfo -*-\n\
00775 @deftypefn {Loadable Function} {[@var{Q1}, @var{R1}] =} qrupdate (@var{Q}, @var{R}, @var{u}, @var{v})\n\
00776 Given a QR@tie{}factorization of a real or complex matrix\n\
00777 @w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and\n\
00778 @var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization\n\
00779 of @w{@var{A} + @var{u}*@var{v}'}, where @var{u} and @var{v} are\n\
00780 column vectors (rank-1 update) or matrices with equal number of columns\n\
00781 (rank-k update).  Notice that the latter case is done as a sequence of rank-1\n\
00782 updates; thus, for k large enough, it will be both faster and more accurate\n\
00783 to recompute the factorization from scratch.\n\
00784 \n\
00785 The QR@tie{}factorization supplied may be either full\n\
00786 (Q is square) or economized (R is square).\n\
00787 \n\
00788 @seealso{qr, qrinsert, qrdelete}\n\
00789 @end deftypefn")
00790 {
00791   octave_idx_type nargin = args.length ();
00792   octave_value_list retval;
00793 
00794   if (nargin != 4)
00795     {
00796       print_usage ();
00797       return retval;
00798     }
00799 
00800   octave_value argq = args(0);
00801   octave_value argr = args(1);
00802   octave_value argu = args(2);
00803   octave_value argv = args(3);
00804 
00805   if (argq.is_numeric_type () && argr.is_numeric_type ()
00806       && argu.is_numeric_type () && argv.is_numeric_type ())
00807     {
00808       if (check_qr_dims (argq, argr, true))
00809         {
00810           if (argq.is_real_type ()
00811               && argr.is_real_type ()
00812               && argu.is_real_type ()
00813               && argv.is_real_type ())
00814             {
00815               // all real case
00816               if (argq.is_single_type ()
00817                   || argr.is_single_type ()
00818                   || argu.is_single_type ()
00819                   || argv.is_single_type ())
00820                 {
00821                   FloatMatrix Q = argq.float_matrix_value ();
00822                   FloatMatrix R = argr.float_matrix_value ();
00823                   FloatMatrix u = argu.float_matrix_value ();
00824                   FloatMatrix v = argv.float_matrix_value ();
00825 
00826                   FloatQR fact (Q, R);
00827                   fact.update (u, v);
00828 
00829                   retval(1) = get_qr_r (fact);
00830                   retval(0) = fact.Q ();
00831                 }
00832               else
00833                 {
00834                   Matrix Q = argq.matrix_value ();
00835                   Matrix R = argr.matrix_value ();
00836                   Matrix u = argu.matrix_value ();
00837                   Matrix v = argv.matrix_value ();
00838 
00839                   QR fact (Q, R);
00840                   fact.update (u, v);
00841 
00842                   retval(1) = get_qr_r (fact);
00843                   retval(0) = fact.Q ();
00844                 }
00845             }
00846           else
00847             {
00848               // complex case
00849               if (argq.is_single_type ()
00850                   || argr.is_single_type ()
00851                   || argu.is_single_type ()
00852                   || argv.is_single_type ())
00853                 {
00854                   FloatComplexMatrix Q = argq.float_complex_matrix_value ();
00855                   FloatComplexMatrix R = argr.float_complex_matrix_value ();
00856                   FloatComplexMatrix u = argu.float_complex_matrix_value ();
00857                   FloatComplexMatrix v = argv.float_complex_matrix_value ();
00858 
00859                   FloatComplexQR fact (Q, R);
00860                   fact.update (u, v);
00861 
00862                   retval(1) = get_qr_r (fact);
00863                   retval(0) = fact.Q ();
00864                 }
00865               else
00866                 {
00867                   ComplexMatrix Q = argq.complex_matrix_value ();
00868                   ComplexMatrix R = argr.complex_matrix_value ();
00869                   ComplexMatrix u = argu.complex_matrix_value ();
00870                   ComplexMatrix v = argv.complex_matrix_value ();
00871 
00872                   ComplexQR fact (Q, R);
00873                   fact.update (u, v);
00874 
00875                   retval(1) = get_qr_r (fact);
00876                   retval(0) = fact.Q ();
00877                 }
00878             }
00879         }
00880       else
00881         error ("qrupdate: Q and R dimensions don't match");
00882     }
00883   else
00884     error ("qrupdate: Q, R, U, and V must be numeric");
00885 
00886   return retval;
00887 }
00888 /*
00889 %!shared A, u, v, Ac, uc, vc
00890 %! A = [0.091364  0.613038  0.999083;
00891 %!      0.594638  0.425302  0.603537;
00892 %!      0.383594  0.291238  0.085574;
00893 %!      0.265712  0.268003  0.238409;
00894 %!      0.669966  0.743851  0.445057 ];
00895 %!
00896 %! u = [0.85082;
00897 %!      0.76426;
00898 %!      0.42883;
00899 %!      0.53010;
00900 %!      0.80683 ];
00901 %!
00902 %! v = [0.98810;
00903 %!      0.24295;
00904 %!      0.43167 ];
00905 %!
00906 %! Ac = [0.620405 + 0.956953i  0.480013 + 0.048806i  0.402627 + 0.338171i;
00907 %!      0.589077 + 0.658457i  0.013205 + 0.279323i  0.229284 + 0.721929i;
00908 %!      0.092758 + 0.345687i  0.928679 + 0.241052i  0.764536 + 0.832406i;
00909 %!      0.912098 + 0.721024i  0.049018 + 0.269452i  0.730029 + 0.796517i;
00910 %!      0.112849 + 0.603871i  0.486352 + 0.142337i  0.355646 + 0.151496i ];
00911 %!
00912 %! uc = [0.20351 + 0.05401i;
00913 %!      0.13141 + 0.43708i;
00914 %!      0.29808 + 0.08789i;
00915 %!      0.69821 + 0.38844i;
00916 %!      0.74871 + 0.25821i ];
00917 %!
00918 %! vc = [0.85839 + 0.29468i;
00919 %!      0.20820 + 0.93090i;
00920 %!      0.86184 + 0.34689i ];
00921 %!
00922 
00923 %!test
00924 %! [Q,R] = qr(A);
00925 %! [Q,R] = qrupdate(Q,R,u,v);
00926 %! assert(norm(vec(Q'*Q - eye(5)),Inf) < 1e1*eps)
00927 %! assert(norm(vec(triu(R)-R),Inf) == 0)
00928 %! assert(norm(vec(Q*R - A - u*v'),Inf) < norm(A)*1e1*eps)
00929 %!
00930 %!test
00931 %! [Q,R] = qr(Ac);
00932 %! [Q,R] = qrupdate(Q,R,uc,vc);
00933 %! assert(norm(vec(Q'*Q - eye(5)),Inf) < 1e1*eps)
00934 %! assert(norm(vec(triu(R)-R),Inf) == 0)
00935 %! assert(norm(vec(Q*R - Ac - uc*vc'),Inf) < norm(Ac)*1e1*eps)
00936 
00937 %!test
00938 %! [Q,R] = qr(single(A));
00939 %! [Q,R] = qrupdate(Q,R,single(u),single(v));
00940 %! assert(norm(vec(Q'*Q - eye(5,'single')),Inf) < 1e1*eps('single'))
00941 %! assert(norm(vec(triu(R)-R),Inf) == 0)
00942 %! assert(norm(vec(Q*R - single(A) - single(u)*single(v)'),Inf) < norm(single(A))*1e1*eps('single'))
00943 %!
00944 %!test
00945 %! [Q,R] = qr(single(Ac));
00946 %! [Q,R] = qrupdate(Q,R,single(uc),single(vc));
00947 %! assert(norm(vec(Q'*Q - eye(5,'single')),Inf) < 1e1*eps('single'))
00948 %! assert(norm(vec(triu(R)-R),Inf) == 0)
00949 %! assert(norm(vec(Q*R - single(Ac) - single(uc)*single(vc)'),Inf) < norm(single(Ac))*1e1*eps('single'))
00950 */
00951 
00952 DEFUN_DLD (qrinsert, args, ,
00953   "-*- texinfo -*-\n\
00954 @deftypefn {Loadable Function} {[@var{Q1}, @var{R1}] =} qrinsert (@var{Q}, @var{R}, @var{j}, @var{x}, @var{orient})\n\
00955 Given a QR@tie{}factorization of a real or complex matrix\n\
00956 @w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and\n\
00957 @var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization of\n\
00958 @w{[A(:,1:j-1) x A(:,j:n)]}, where @var{u} is a column vector to be\n\
00959 inserted into @var{A} (if @var{orient} is @code{\"col\"}), or the\n\
00960 QR@tie{}factorization of @w{[A(1:j-1,:);x;A(:,j:n)]}, where @var{x}\n\
00961 is a row vector to be inserted into @var{A} (if @var{orient} is\n\
00962 @code{\"row\"}).\n\
00963 \n\
00964 The default value of @var{orient} is @code{\"col\"}.\n\
00965 If @var{orient} is @code{\"col\"},\n\
00966 @var{u} may be a matrix and @var{j} an index vector\n\
00967 resulting in the QR@tie{}factorization of a matrix @var{B} such that\n\
00968 @w{B(:,@var{j})} gives @var{u} and @w{B(:,@var{j}) = []} gives @var{A}.\n\
00969 Notice that the latter case is done as a sequence of k insertions;\n\
00970 thus, for k large enough, it will be both faster and more accurate to\n\
00971 recompute the factorization from scratch.\n\
00972 \n\
00973 If @var{orient} is @code{\"col\"},\n\
00974 the QR@tie{}factorization supplied may be either full\n\
00975 (Q is square) or economized (R is square).\n\
00976 \n\
00977 If @var{orient} is @code{\"row\"}, full factorization is needed.\n\
00978 @seealso{qr, qrupdate, qrdelete}\n\
00979 @end deftypefn")
00980 {
00981   octave_idx_type nargin = args.length ();
00982   octave_value_list retval;
00983 
00984   if (nargin < 4 || nargin > 5)
00985     {
00986       print_usage ();
00987       return retval;
00988     }
00989 
00990   octave_value argq = args(0);
00991   octave_value argr = args(1);
00992   octave_value argj = args(2);
00993   octave_value argx = args(3);
00994 
00995   if (argq.is_numeric_type () && argr.is_numeric_type ()
00996       && argx.is_numeric_type ()
00997       && (nargin < 5 || args(4).is_string ()))
00998     {
00999       std::string orient = (nargin < 5) ? "col" : args(4).string_value ();
01000 
01001       bool col = orient == "col";
01002 
01003       if (col || orient == "row")
01004         if (check_qr_dims (argq, argr, col)
01005             && (col || argx.rows () == 1))
01006           {
01007             if (check_index (argj, col))
01008               {
01009                 MArray<octave_idx_type> j
01010                   = argj.octave_idx_type_vector_value ();
01011 
01012                 octave_idx_type one = 1;
01013 
01014                 if (argq.is_real_type ()
01015                     && argr.is_real_type ()
01016                     && argx.is_real_type ())
01017                   {
01018                     // real case
01019                     if (argq.is_single_type ()
01020                         || argr.is_single_type ()
01021                         || argx.is_single_type ())
01022                       {
01023                         FloatMatrix Q = argq.float_matrix_value ();
01024                         FloatMatrix R = argr.float_matrix_value ();
01025                         FloatMatrix x = argx.float_matrix_value ();
01026 
01027                         FloatQR fact (Q, R);
01028 
01029                         if (col)
01030                           fact.insert_col (x, j-one);
01031                         else
01032                           fact.insert_row (x.row (0), j(0)-one);
01033 
01034                         retval(1) = get_qr_r (fact);
01035                         retval(0) = fact.Q ();
01036 
01037                       }
01038                     else
01039                       {
01040                         Matrix Q = argq.matrix_value ();
01041                         Matrix R = argr.matrix_value ();
01042                         Matrix x = argx.matrix_value ();
01043 
01044                         QR fact (Q, R);
01045 
01046                         if (col)
01047                           fact.insert_col (x, j-one);
01048                         else
01049                           fact.insert_row (x.row (0), j(0)-one);
01050 
01051                         retval(1) = get_qr_r (fact);
01052                         retval(0) = fact.Q ();
01053 
01054                       }
01055                   }
01056                 else
01057                   {
01058                     // complex case
01059                     if (argq.is_single_type ()
01060                         || argr.is_single_type ()
01061                         || argx.is_single_type ())
01062                       {
01063                         FloatComplexMatrix Q = argq.float_complex_matrix_value ();
01064                         FloatComplexMatrix R = argr.float_complex_matrix_value ();
01065                         FloatComplexMatrix x = argx.float_complex_matrix_value ();
01066 
01067                         FloatComplexQR fact (Q, R);
01068 
01069                         if (col)
01070                           fact.insert_col (x, j-one);
01071                         else
01072                           fact.insert_row (x.row (0), j(0)-one);
01073 
01074                         retval(1) = get_qr_r (fact);
01075                         retval(0) = fact.Q ();
01076                       }
01077                     else
01078                       {
01079                         ComplexMatrix Q = argq.complex_matrix_value ();
01080                         ComplexMatrix R = argr.complex_matrix_value ();
01081                         ComplexMatrix x = argx.complex_matrix_value ();
01082 
01083                         ComplexQR fact (Q, R);
01084 
01085                         if (col)
01086                           fact.insert_col (x, j-one);
01087                         else
01088                           fact.insert_row (x.row (0), j(0)-one);
01089 
01090                         retval(1) = get_qr_r (fact);
01091                         retval(0) = fact.Q ();
01092                       }
01093                   }
01094 
01095               }
01096             else
01097               error ("qrinsert: invalid index J");
01098           }
01099         else
01100           error ("qrinsert: dimension mismatch");
01101 
01102       else
01103         error ("qrinsert: ORIENT must be \"col\" or \"row\"");
01104     }
01105   else
01106     print_usage ();
01107 
01108   return retval;
01109 }
01110 
01111 /*
01112 %!test
01113 %! [Q,R] = qr(A);
01114 %! [Q,R] = qrinsert(Q,R,3,u);
01115 %! assert(norm(vec(Q'*Q - eye(5)),Inf) < 1e1*eps)
01116 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01117 %! assert(norm(vec(Q*R - [A(:,1:2) u A(:,3)]),Inf) < norm(A)*1e1*eps)
01118 %!test
01119 %! [Q,R] = qr(Ac);
01120 %! [Q,R] = qrinsert(Q,R,3,uc);
01121 %! assert(norm(vec(Q'*Q - eye(5)),Inf) < 1e1*eps)
01122 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01123 %! assert(norm(vec(Q*R - [Ac(:,1:2) uc Ac(:,3)]),Inf) < norm(Ac)*1e1*eps)
01124 %!test
01125 %! x = [0.85082  0.76426  0.42883 ];
01126 %!
01127 %! [Q,R] = qr(A);
01128 %! [Q,R] = qrinsert(Q,R,3,x,'row');
01129 %! assert(norm(vec(Q'*Q - eye(6)),Inf) < 1e1*eps)
01130 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01131 %! assert(norm(vec(Q*R - [A(1:2,:);x;A(3:5,:)]),Inf) < norm(A)*1e1*eps)
01132 %!test
01133 %! x = [0.20351 + 0.05401i  0.13141 + 0.43708i  0.29808 + 0.08789i ];
01134 %!
01135 %! [Q,R] = qr(Ac);
01136 %! [Q,R] = qrinsert(Q,R,3,x,'row');
01137 %! assert(norm(vec(Q'*Q - eye(6)),Inf) < 1e1*eps)
01138 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01139 %! assert(norm(vec(Q*R - [Ac(1:2,:);x;Ac(3:5,:)]),Inf) < norm(Ac)*1e1*eps)
01140 
01141 %!test
01142 %! [Q,R] = qr(single(A));
01143 %! [Q,R] = qrinsert(Q,R,3,single(u));
01144 %! assert(norm(vec(Q'*Q - eye(5,'single')),Inf) < 1e1*eps('single'))
01145 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01146 %! assert(norm(vec(Q*R - single([A(:,1:2) u A(:,3)])),Inf) < norm(single(A))*1e1*eps('single'))
01147 %!test
01148 %! [Q,R] = qr(single(Ac));
01149 %! [Q,R] = qrinsert(Q,R,3,single(uc));
01150 %! assert(norm(vec(Q'*Q - eye(5,'single')),Inf) < 1e1*eps('single'))
01151 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01152 %! assert(norm(vec(Q*R - single([Ac(:,1:2) uc Ac(:,3)])),Inf) < norm(single(Ac))*1e1*eps('single'))
01153 %!test
01154 %! x = single([0.85082  0.76426  0.42883 ]);
01155 %!
01156 %! [Q,R] = qr(single(A));
01157 %! [Q,R] = qrinsert(Q,R,3,x,'row');
01158 %! assert(norm(vec(Q'*Q - eye(6,'single')),Inf) < 1e1*eps('single'))
01159 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01160 %! assert(norm(vec(Q*R - single([A(1:2,:);x;A(3:5,:)])),Inf) < norm(single(A))*1e1*eps('single'))
01161 %!test
01162 %! x = single([0.20351 + 0.05401i  0.13141 + 0.43708i  0.29808 + 0.08789i ]);
01163 %!
01164 %! [Q,R] = qr(single(Ac));
01165 %! [Q,R] = qrinsert(Q,R,3,x,'row');
01166 %! assert(norm(vec(Q'*Q - eye(6,'single')),Inf) < 1e1*eps('single'))
01167 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01168 %! assert(norm(vec(Q*R - single([Ac(1:2,:);x;Ac(3:5,:)])),Inf) < norm(single(Ac))*1e1*eps('single'))
01169 */
01170 
01171 DEFUN_DLD (qrdelete, args, ,
01172   "-*- texinfo -*-\n\
01173 @deftypefn {Loadable Function} {[@var{Q1}, @var{R1}] =} qrdelete (@var{Q}, @var{R}, @var{j}, @var{orient})\n\
01174 Given a QR@tie{}factorization of a real or complex matrix\n\
01175 @w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and\n\
01176 @var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization of\n\
01177 @w{[A(:,1:j-1) A(:,j+1:n)]}, i.e., @var{A} with one column deleted\n\
01178 (if @var{orient} is \"col\"), or the QR@tie{}factorization of\n\
01179 @w{[A(1:j-1,:);A(j+1:n,:)]}, i.e., @var{A} with one row deleted (if\n   \
01180 @var{orient} is \"row\").\n\
01181 \n\
01182 The default value of @var{orient} is \"col\".\n\
01183 \n\
01184 If @var{orient} is @code{\"col\"},\n\
01185 @var{j} may be an index vector\n\
01186 resulting in the QR@tie{}factorization of a matrix @var{B} such that\n\
01187 @w{A(:,@var{j}) = []} gives @var{B}.\n\
01188 Notice that the latter case is done as a sequence of k deletions;\n\
01189 thus, for k large enough, it will be both faster and more accurate to\n\
01190 recompute the factorization from scratch.\n\
01191 \n\
01192 If @var{orient} is @code{\"col\"},\n\
01193 the QR@tie{}factorization supplied may be either full\n\
01194 (Q is square) or economized (R is square).\n\
01195 \n\
01196 If @var{orient} is @code{\"row\"}, full factorization is needed.\n\
01197 @seealso{qr, qrinsert, qrupdate}\n\
01198 @end deftypefn")
01199 {
01200   octave_idx_type nargin = args.length ();
01201   octave_value_list retval;
01202 
01203   if (nargin < 3 || nargin > 4)
01204     {
01205       print_usage ();
01206       return retval;
01207     }
01208 
01209   octave_value argq = args(0);
01210   octave_value argr = args(1);
01211   octave_value argj = args(2);
01212 
01213   if (argq.is_numeric_type () && argr.is_numeric_type ()
01214       && (nargin < 4 || args(3).is_string ()))
01215     {
01216       std::string orient = (nargin < 4) ? "col" : args(3).string_value ();
01217 
01218       bool col = orient == "col";
01219 
01220       if (col || orient == "row")
01221         if (check_qr_dims (argq, argr, col))
01222           {
01223             if (check_index (argj, col))
01224               {
01225                 MArray<octave_idx_type> j
01226                   = argj.octave_idx_type_vector_value ();
01227 
01228                 octave_idx_type one = 1;
01229 
01230                 if (argq.is_real_type ()
01231                     && argr.is_real_type ())
01232                   {
01233                     // real case
01234                     if (argq.is_single_type ()
01235                         || argr.is_single_type ())
01236                       {
01237                         FloatMatrix Q = argq.float_matrix_value ();
01238                         FloatMatrix R = argr.float_matrix_value ();
01239 
01240                         FloatQR fact (Q, R);
01241 
01242                         if (col)
01243                           fact.delete_col (j-one);
01244                         else
01245                           fact.delete_row (j(0)-one);
01246 
01247                         retval(1) = get_qr_r (fact);
01248                         retval(0) = fact.Q ();
01249                       }
01250                     else
01251                       {
01252                         Matrix Q = argq.matrix_value ();
01253                         Matrix R = argr.matrix_value ();
01254 
01255                         QR fact (Q, R);
01256 
01257                         if (col)
01258                           fact.delete_col (j-one);
01259                         else
01260                           fact.delete_row (j(0)-one);
01261 
01262                         retval(1) = get_qr_r (fact);
01263                         retval(0) = fact.Q ();
01264                       }
01265                   }
01266                 else
01267                   {
01268                     // complex case
01269                     if (argq.is_single_type ()
01270                         || argr.is_single_type ())
01271                       {
01272                         FloatComplexMatrix Q = argq.float_complex_matrix_value ();
01273                         FloatComplexMatrix R = argr.float_complex_matrix_value ();
01274 
01275                         FloatComplexQR fact (Q, R);
01276 
01277                         if (col)
01278                           fact.delete_col (j-one);
01279                         else
01280                           fact.delete_row (j(0)-one);
01281 
01282                         retval(1) = get_qr_r (fact);
01283                         retval(0) = fact.Q ();
01284                       }
01285                     else
01286                       {
01287                         ComplexMatrix Q = argq.complex_matrix_value ();
01288                         ComplexMatrix R = argr.complex_matrix_value ();
01289 
01290                         ComplexQR fact (Q, R);
01291 
01292                         if (col)
01293                           fact.delete_col (j-one);
01294                         else
01295                           fact.delete_row (j(0)-one);
01296 
01297                         retval(1) = get_qr_r (fact);
01298                         retval(0) = fact.Q ();
01299                       }
01300                   }
01301               }
01302             else
01303               error ("qrdelete: invalid index J");
01304           }
01305         else
01306           error ("qrdelete: dimension mismatch");
01307 
01308       else
01309         error ("qrdelete: ORIENT must be \"col\" or \"row\"");
01310     }
01311   else
01312     print_usage ();
01313 
01314   return retval;
01315 }
01316 
01317 /*
01318 %!test
01319 %! AA = [0.091364  0.613038  0.027504  0.999083;
01320 %!       0.594638  0.425302  0.562834  0.603537;
01321 %!       0.383594  0.291238  0.742073  0.085574;
01322 %!       0.265712  0.268003  0.783553  0.238409;
01323 %!       0.669966  0.743851  0.457255  0.445057 ];
01324 %!
01325 %! [Q,R] = qr(AA);
01326 %! [Q,R] = qrdelete(Q,R,3);
01327 %! assert(norm(vec(Q'*Q - eye(5)),Inf) < 16*eps)
01328 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01329 %! assert(norm(vec(Q*R - [AA(:,1:2) AA(:,4)]),Inf) < norm(AA)*1e1*eps)
01330 %!
01331 %!test
01332 %! AA = [0.364554 + 0.993117i  0.669818 + 0.510234i  0.426568 + 0.041337i  0.847051 + 0.233291i;
01333 %!       0.049600 + 0.242783i  0.448946 + 0.484022i  0.141155 + 0.074420i  0.446746 + 0.392706i;
01334 %!       0.581922 + 0.657416i  0.581460 + 0.030016i  0.219909 + 0.447288i  0.201144 + 0.069132i;
01335 %!       0.694986 + 0.000571i  0.682327 + 0.841712i  0.807537 + 0.166086i  0.192767 + 0.358098i;
01336 %!       0.945002 + 0.066788i  0.350492 + 0.642638i  0.579629 + 0.048102i  0.600170 + 0.636938i ] * I;
01337 %!
01338 %! [Q,R] = qr(AA);
01339 %! [Q,R] = qrdelete(Q,R,3);
01340 %! assert(norm(vec(Q'*Q - eye(5)),Inf) < 16*eps)
01341 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01342 %! assert(norm(vec(Q*R - [AA(:,1:2) AA(:,4)]),Inf) < norm(AA)*1e1*eps)
01343 %!
01344 %!test
01345 %! AA = [0.091364  0.613038  0.027504  0.999083;
01346 %!       0.594638  0.425302  0.562834  0.603537;
01347 %!       0.383594  0.291238  0.742073  0.085574;
01348 %!       0.265712  0.268003  0.783553  0.238409;
01349 %!       0.669966  0.743851  0.457255  0.445057 ];
01350 %!
01351 %! [Q,R] = qr(AA);
01352 %! [Q,R] = qrdelete(Q,R,3,'row');
01353 %! assert(norm(vec(Q'*Q - eye(4)),Inf) < 1e1*eps)
01354 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01355 %! assert(norm(vec(Q*R - [AA(1:2,:);AA(4:5,:)]),Inf) < norm(AA)*1e1*eps)
01356 %!
01357 %!test
01358 %! AA = [0.364554 + 0.993117i  0.669818 + 0.510234i  0.426568 + 0.041337i  0.847051 + 0.233291i;
01359 %!       0.049600 + 0.242783i  0.448946 + 0.484022i  0.141155 + 0.074420i  0.446746 + 0.392706i;
01360 %!       0.581922 + 0.657416i  0.581460 + 0.030016i  0.219909 + 0.447288i  0.201144 + 0.069132i;
01361 %!       0.694986 + 0.000571i  0.682327 + 0.841712i  0.807537 + 0.166086i  0.192767 + 0.358098i;
01362 %!       0.945002 + 0.066788i  0.350492 + 0.642638i  0.579629 + 0.048102i  0.600170 + 0.636938i ] * I;
01363 %!
01364 %! [Q,R] = qr(AA);
01365 %! [Q,R] = qrdelete(Q,R,3,'row');
01366 %! assert(norm(vec(Q'*Q - eye(4)),Inf) < 1e1*eps)
01367 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01368 %! assert(norm(vec(Q*R - [AA(1:2,:);AA(4:5,:)]),Inf) < norm(AA)*1e1*eps)
01369 
01370 %!test
01371 %! AA = single([0.091364  0.613038  0.027504  0.999083;
01372 %!              0.594638  0.425302  0.562834  0.603537;
01373 %!              0.383594  0.291238  0.742073  0.085574;
01374 %!              0.265712  0.268003  0.783553  0.238409;
01375 %!              0.669966  0.743851  0.457255  0.445057 ]);
01376 %!
01377 %! [Q,R] = qr(AA);
01378 %! [Q,R] = qrdelete(Q,R,3);
01379 %! assert(norm(vec(Q'*Q - eye(5,'single')),Inf) < 1e1*eps('single'))
01380 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01381 %! assert(norm(vec(Q*R - [AA(:,1:2) AA(:,4)]),Inf) < norm(AA)*1e1*eps('single'))
01382 %!
01383 %!test
01384 %! AA = single([0.364554 + 0.993117i  0.669818 + 0.510234i  0.426568 + 0.041337i  0.847051 + 0.233291i;
01385 %!              0.049600 + 0.242783i  0.448946 + 0.484022i  0.141155 + 0.074420i  0.446746 + 0.392706i;
01386 %!              0.581922 + 0.657416i  0.581460 + 0.030016i  0.219909 + 0.447288i  0.201144 + 0.069132i;
01387 %!              0.694986 + 0.000571i  0.682327 + 0.841712i  0.807537 + 0.166086i  0.192767 + 0.358098i;
01388 %!              0.945002 + 0.066788i  0.350492 + 0.642638i  0.579629 + 0.048102i  0.600170 + 0.636938i ]) * I;
01389 %!
01390 %! [Q,R] = qr(AA);
01391 %! [Q,R] = qrdelete(Q,R,3);
01392 %! assert(norm(vec(Q'*Q - eye(5,'single')),Inf) < 1e1*eps('single'))
01393 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01394 %! assert(norm(vec(Q*R - [AA(:,1:2) AA(:,4)]),Inf) < norm(AA)*1e1*eps('single'))
01395 %!
01396 %!test
01397 %! AA = single([0.091364  0.613038  0.027504  0.999083;
01398 %!              0.594638  0.425302  0.562834  0.603537;
01399 %!              0.383594  0.291238  0.742073  0.085574;
01400 %!              0.265712  0.268003  0.783553  0.238409;
01401 %!              0.669966  0.743851  0.457255  0.445057 ]);
01402 %!
01403 %! [Q,R] = qr(AA);
01404 %! [Q,R] = qrdelete(Q,R,3,'row');
01405 %! assert(norm(vec(Q'*Q - eye(4,'single')),Inf) < 1.5e1*eps('single'))
01406 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01407 %! assert(norm(vec(Q*R - [AA(1:2,:);AA(4:5,:)]),Inf) < norm(AA)*1e1*eps('single'))
01408 %!testif HAVE_QRUPDATE
01409 %! # Same test as above but with more precicision
01410 %! AA = single([0.091364  0.613038  0.027504  0.999083;
01411 %!              0.594638  0.425302  0.562834  0.603537;
01412 %!              0.383594  0.291238  0.742073  0.085574;
01413 %!              0.265712  0.268003  0.783553  0.238409;
01414 %!              0.669966  0.743851  0.457255  0.445057 ]);
01415 %!
01416 %! [Q,R] = qr(AA);
01417 %! [Q,R] = qrdelete(Q,R,3,'row');
01418 %! assert(norm(vec(Q'*Q - eye(4,'single')),Inf) < 1e1*eps('single'))
01419 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01420 %! assert(norm(vec(Q*R - [AA(1:2,:);AA(4:5,:)]),Inf) < norm(AA)*1e1*eps('single'))
01421 %!
01422 %!test
01423 %! AA = single([0.364554 + 0.993117i  0.669818 + 0.510234i  0.426568 + 0.041337i  0.847051 + 0.233291i;
01424 %!              0.049600 + 0.242783i  0.448946 + 0.484022i  0.141155 + 0.074420i  0.446746 + 0.392706i;
01425 %!              0.581922 + 0.657416i  0.581460 + 0.030016i  0.219909 + 0.447288i  0.201144 + 0.069132i;
01426 %!              0.694986 + 0.000571i  0.682327 + 0.841712i  0.807537 + 0.166086i  0.192767 + 0.358098i;
01427 %!              0.945002 + 0.066788i  0.350492 + 0.642638i  0.579629 + 0.048102i  0.600170 + 0.636938i ]) * I;
01428 %!
01429 %! [Q,R] = qr(AA);
01430 %! [Q,R] = qrdelete(Q,R,3,'row');
01431 %! assert(norm(vec(Q'*Q - eye(4,'single')),Inf) < 1e1*eps('single'))
01432 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01433 %! assert(norm(vec(Q*R - [AA(1:2,:);AA(4:5,:)]),Inf) < norm(AA)*1e1*eps('single'))
01434 */
01435 
01436 DEFUN_DLD (qrshift, args, ,
01437   "-*- texinfo -*-\n\
01438 @deftypefn {Loadable Function} {[@var{Q1}, @var{R1}] =} qrshift (@var{Q}, @var{R}, @var{i}, @var{j})\n\
01439 Given a QR@tie{}factorization of a real or complex matrix\n\
01440 @w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and\n\
01441 @var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization\n\
01442 of @w{@var{A}(:,p)}, where @w{p} is the permutation @*\n\
01443 @code{p = [1:i-1, shift(i:j, 1), j+1:n]} if @w{@var{i} < @var{j}} @*\n\
01444  or @*\n\
01445 @code{p = [1:j-1, shift(j:i,-1), i+1:n]} if @w{@var{j} < @var{i}}.  @*\n\
01446 \n\
01447 @seealso{qr, qrinsert, qrdelete}\n\
01448 @end deftypefn")
01449 {
01450   octave_idx_type nargin = args.length ();
01451   octave_value_list retval;
01452 
01453   if (nargin != 4)
01454     {
01455       print_usage ();
01456       return retval;
01457     }
01458 
01459   octave_value argq = args(0);
01460   octave_value argr = args(1);
01461   octave_value argi = args(2);
01462   octave_value argj = args(3);
01463 
01464   if (argq.is_numeric_type () && argr.is_numeric_type ())
01465     {
01466       if (check_qr_dims (argq, argr, true))
01467         {
01468           if (check_index (argi) && check_index (argj))
01469             {
01470               octave_idx_type i = argi.int_value ();
01471               octave_idx_type j = argj.int_value ();
01472 
01473               if (argq.is_real_type ()
01474                   && argr.is_real_type ())
01475                 {
01476                   // all real case
01477                   if (argq.is_single_type ()
01478                       && argr.is_single_type ())
01479                     {
01480                       FloatMatrix Q = argq.float_matrix_value ();
01481                       FloatMatrix R = argr.float_matrix_value ();
01482 
01483                       FloatQR fact (Q, R);
01484                       fact.shift_cols (i-1, j-1);
01485 
01486                       retval(1) = get_qr_r (fact);
01487                       retval(0) = fact.Q ();
01488                     }
01489                   else
01490                     {
01491                       Matrix Q = argq.matrix_value ();
01492                       Matrix R = argr.matrix_value ();
01493 
01494                       QR fact (Q, R);
01495                       fact.shift_cols (i-1, j-1);
01496 
01497                       retval(1) = get_qr_r (fact);
01498                       retval(0) = fact.Q ();
01499                     }
01500                 }
01501               else
01502                 {
01503                   // complex case
01504                   if (argq.is_single_type ()
01505                       && argr.is_single_type ())
01506                     {
01507                       FloatComplexMatrix Q = argq.float_complex_matrix_value ();
01508                       FloatComplexMatrix R = argr.float_complex_matrix_value ();
01509 
01510                       FloatComplexQR fact (Q, R);
01511                       fact.shift_cols (i-1, j-1);
01512 
01513                       retval(1) = get_qr_r (fact);
01514                       retval(0) = fact.Q ();
01515                     }
01516                   else
01517                     {
01518                       ComplexMatrix Q = argq.complex_matrix_value ();
01519                       ComplexMatrix R = argr.complex_matrix_value ();
01520 
01521                       ComplexQR fact (Q, R);
01522                       fact.shift_cols (i-1, j-1);
01523 
01524                       retval(1) = get_qr_r (fact);
01525                       retval(0) = fact.Q ();
01526                     }
01527                 }
01528             }
01529           else
01530             error ("qrshift: invalid index I or J");
01531         }
01532       else
01533         error ("qrshift: dimensions mismatch");
01534     }
01535   else
01536     error ("qrshift: Q and R must be numeric");
01537 
01538   return retval;
01539 }
01540 /*
01541 %!test
01542 %! AA = A.';
01543 %! i = 2; j = 4; p = [1:i-1, shift(i:j,-1), j+1:5];
01544 %!
01545 %! [Q,R] = qr(AA);
01546 %! [Q,R] = qrshift(Q,R,i,j);
01547 %! assert(norm(vec(Q'*Q - eye(3)),Inf) < 1e1*eps)
01548 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01549 %! assert(norm(vec(Q*R - AA(:,p)),Inf) < norm(AA)*1e1*eps)
01550 %!
01551 %! j = 2; i = 4; p = [1:j-1, shift(j:i,+1), i+1:5];
01552 %!
01553 %! [Q,R] = qr(AA);
01554 %! [Q,R] = qrshift(Q,R,i,j);
01555 %! assert(norm(vec(Q'*Q - eye(3)),Inf) < 1e1*eps)
01556 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01557 %! assert(norm(vec(Q*R - AA(:,p)),Inf) < norm(AA)*1e1*eps)
01558 %!
01559 %!test
01560 %! AA = Ac.';
01561 %! i = 2; j = 4; p = [1:i-1, shift(i:j,-1), j+1:5];
01562 %!
01563 %! [Q,R] = qr(AA);
01564 %! [Q,R] = qrshift(Q,R,i,j);
01565 %! assert(norm(vec(Q'*Q - eye(3)),Inf) < 1e1*eps)
01566 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01567 %! assert(norm(vec(Q*R - AA(:,p)),Inf) < norm(AA)*1e1*eps)
01568 %!
01569 %! j = 2; i = 4; p = [1:j-1, shift(j:i,+1), i+1:5];
01570 %!
01571 %! [Q,R] = qr(AA);
01572 %! [Q,R] = qrshift(Q,R,i,j);
01573 %! assert(norm(vec(Q'*Q - eye(3)),Inf) < 1e1*eps)
01574 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01575 %! assert(norm(vec(Q*R - AA(:,p)),Inf) < norm(AA)*1e1*eps)
01576 
01577 
01578 %!test
01579 %! AA = single (A).';
01580 %! i = 2; j = 4; p = [1:i-1, shift(i:j,-1), j+1:5];
01581 %!
01582 %! [Q,R] = qr(AA);
01583 %! [Q,R] = qrshift(Q,R,i,j);
01584 %! assert(norm(vec(Q'*Q - eye(3,'single')),Inf) < 1e1*eps('single'))
01585 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01586 %! assert(norm(vec(Q*R - AA(:,p)),Inf) < norm(AA)*1e1*eps('single'))
01587 %!
01588 %! j = 2; i = 4; p = [1:j-1, shift(j:i,+1), i+1:5];
01589 %!
01590 %! [Q,R] = qr(AA);
01591 %! [Q,R] = qrshift(Q,R,i,j);
01592 %! assert(norm(vec(Q'*Q - eye(3,'single')),Inf) < 1e1*eps('single'))
01593 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01594 %! assert(norm(vec(Q*R - AA(:,p)),Inf) < norm(AA)*1e1*eps('single'))
01595 %!
01596 %!test
01597 %! AA = single(Ac).';
01598 %! i = 2; j = 4; p = [1:i-1, shift(i:j,-1), j+1:5];
01599 %!
01600 %! [Q,R] = qr(AA);
01601 %! [Q,R] = qrshift(Q,R,i,j);
01602 %! assert(norm(vec(Q'*Q - eye(3,'single')),Inf) < 1e1*eps('single'))
01603 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01604 %! assert(norm(vec(Q*R - AA(:,p)),Inf) < norm(AA)*1e1*eps('single'))
01605 %!
01606 %! j = 2; i = 4; p = [1:j-1, shift(j:i,+1), i+1:5];
01607 %!
01608 %! [Q,R] = qr(AA);
01609 %! [Q,R] = qrshift(Q,R,i,j);
01610 %! assert(norm(vec(Q'*Q - eye(3,'single')),Inf) < 1e1*eps('single'))
01611 %! assert(norm(vec(triu(R)-R),Inf) == 0)
01612 %! assert(norm(vec(Q*R - AA(:,p)),Inf) < norm(AA)*1e1*eps('single'))
01613 */
 All Classes Files Functions Variables Typedefs Enumerations Enumerator Friends Defines