00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
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
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
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
00461
00462
00463
00464
00465
00466
00467
00468
00469
00470
00471
00472
00473
00474
00475
00476
00477
00478
00479
00480
00481
00482
00483
00484
00485
00486
00487
00488
00489
00490
00491
00492
00493
00494
00495
00496
00497
00498
00499
00500
00501
00502
00503
00504
00505
00506
00507
00508
00509
00510
00511
00512
00513
00514
00515
00516
00517
00518
00519
00520
00521
00522
00523
00524
00525
00526
00527
00528
00529
00530
00531
00532
00533
00534
00535
00536
00537
00538
00539
00540
00541
00542
00543
00544
00545
00546
00547
00548
00549
00550
00551
00552
00553
00554
00555
00556
00557
00558
00559
00560
00561
00562
00563
00564
00565
00566
00567
00568
00569
00570
00571
00572
00573
00574
00575
00576
00577
00578
00579
00580
00581
00582
00583
00584
00585
00586
00587
00588
00589
00590
00591
00592
00593
00594
00595
00596
00597
00598
00599
00600
00601
00602
00603
00604
00605
00606
00607
00608
00609
00610
00611
00612
00613
00614
00615
00616
00617
00618
00619
00620
00621
00622
00623
00624
00625
00626
00627
00628
00629
00630
00631
00632
00633
00634
00635
00636
00637
00638
00639
00640
00641
00642
00643
00644
00645
00646
00647
00648
00649
00650
00651
00652
00653
00654
00655
00656
00657
00658
00659
00660
00661
00662
00663
00664
00665
00666
00667
00668
00669
00670
00671
00672
00673
00674
00675
00676
00677
00678
00679
00680
00681
00682
00683
00684
00685
00686
00687
00688
00689
00690
00691
00692
00693
00694
00695
00696
00697
00698
00699
00700
00701
00702
00703
00704
00705
00706
00707
00708
00709
00710
00711
00712
00713
00714
00715
00716
00717
00718
00719
00720
00721
00722
00723
00724
00725
00726
00727
00728
00729
00730
00731
00732
00733
00734
00735
00736
00737
00738
00739
00740
00741
00742
00743
00744
00745
00746
00747
00748
00749
00750
00751
00752
00753
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
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
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
00890
00891
00892
00893
00894
00895
00896
00897
00898
00899
00900
00901
00902
00903
00904
00905
00906
00907
00908
00909
00910
00911
00912
00913
00914
00915
00916
00917
00918
00919
00920
00921
00922
00923
00924
00925
00926
00927
00928
00929
00930
00931
00932
00933
00934
00935
00936
00937
00938
00939
00940
00941
00942
00943
00944
00945
00946
00947
00948
00949
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
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
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
01113
01114
01115
01116
01117
01118
01119
01120
01121
01122
01123
01124
01125
01126
01127
01128
01129
01130
01131
01132
01133
01134
01135
01136
01137
01138
01139
01140
01141
01142
01143
01144
01145
01146
01147
01148
01149
01150
01151
01152
01153
01154
01155
01156
01157
01158
01159
01160
01161
01162
01163
01164
01165
01166
01167
01168
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
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
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
01319
01320
01321
01322
01323
01324
01325
01326
01327
01328
01329
01330
01331
01332
01333
01334
01335
01336
01337
01338
01339
01340
01341
01342
01343
01344
01345
01346
01347
01348
01349
01350
01351
01352
01353
01354
01355
01356
01357
01358
01359
01360
01361
01362
01363
01364
01365
01366
01367
01368
01369
01370
01371
01372
01373
01374
01375
01376
01377
01378
01379
01380
01381
01382
01383
01384
01385
01386
01387
01388
01389
01390
01391
01392
01393
01394
01395
01396
01397
01398
01399
01400
01401
01402
01403
01404
01405
01406
01407
01408
01409
01410
01411
01412
01413
01414
01415
01416
01417
01418
01419
01420
01421
01422
01423
01424
01425
01426
01427
01428
01429
01430
01431
01432
01433
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
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
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
01542
01543
01544
01545
01546
01547
01548
01549
01550
01551
01552
01553
01554
01555
01556
01557
01558
01559
01560
01561
01562
01563
01564
01565
01566
01567
01568
01569
01570
01571
01572
01573
01574
01575
01576
01577
01578
01579
01580
01581
01582
01583
01584
01585
01586
01587
01588
01589
01590
01591
01592
01593
01594
01595
01596
01597
01598
01599
01600
01601
01602
01603
01604
01605
01606
01607
01608
01609
01610
01611
01612
01613