00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023 #ifdef HAVE_CONFIG_H
00024 #include <config.h>
00025 #endif
00026
00027 #include "CmplxLU.h"
00028 #include "dbleLU.h"
00029 #include "fCmplxLU.h"
00030 #include "floatLU.h"
00031 #include "SparseCmplxLU.h"
00032 #include "SparsedbleLU.h"
00033
00034 #include "defun-dld.h"
00035 #include "error.h"
00036 #include "gripes.h"
00037 #include "oct-obj.h"
00038 #include "utils.h"
00039 #include "ov-re-sparse.h"
00040 #include "ov-cx-sparse.h"
00041
00042 template <class MT>
00043 static octave_value
00044 get_lu_l (const base_lu<MT>& fact)
00045 {
00046 MT L = fact.L ();
00047 if (L.is_square ())
00048 return octave_value (L, MatrixType (MatrixType::Lower));
00049 else
00050 return L;
00051 }
00052
00053 template <class MT>
00054 static octave_value
00055 get_lu_u (const base_lu<MT>& fact)
00056 {
00057 MT U = fact.U ();
00058 if (U.is_square () && fact.regular ())
00059 return octave_value (U, MatrixType (MatrixType::Upper));
00060 else
00061 return U;
00062 }
00063
00064 DEFUN_DLD (lu, args, nargout,
00065 "-*- texinfo -*-\n\
00066 @deftypefn {Loadable Function} {[@var{L}, @var{U}] =} lu (@var{A})\n\
00067 @deftypefnx {Loadable Function} {[@var{L}, @var{U}, @var{P}] =} lu (@var{A})\n\
00068 @deftypefnx {Loadable Function} {[@var{L}, @var{U}, @var{P}, @var{Q}] =} lu (@var{S})\n\
00069 @deftypefnx {Loadable Function} {[@var{L}, @var{U}, @var{P}, @var{Q}, @var{R}] =} lu (@var{S})\n\
00070 @deftypefnx {Loadable Function} {[@dots{}] =} lu (@var{S}, @var{thres})\n\
00071 @deftypefnx {Loadable Function} {@var{y} =} lu (@dots{})\n\
00072 @deftypefnx {Loadable Function} {[@dots{}] =} lu (@dots{}, 'vector')\n\
00073 @cindex LU decomposition\n\
00074 Compute the LU@tie{}decomposition of @var{A}. If @var{A} is full\n\
00075 subroutines from\n\
00076 @sc{lapack} are used and if @var{A} is sparse then @sc{umfpack} is used. The\n\
00077 result is returned in a permuted form, according to the optional return\n\
00078 value @var{P}. For example, given the matrix @code{a = [1, 2; 3, 4]},\n\
00079 \n\
00080 @example\n\
00081 [l, u, p] = lu (@var{a})\n\
00082 @end example\n\
00083 \n\
00084 @noindent\n\
00085 returns\n\
00086 \n\
00087 @example\n\
00088 @group\n\
00089 l =\n\
00090 \n\
00091 1.00000 0.00000\n\
00092 0.33333 1.00000\n\
00093 \n\
00094 u =\n\
00095 \n\
00096 3.00000 4.00000\n\
00097 0.00000 0.66667\n\
00098 \n\
00099 p =\n\
00100 \n\
00101 0 1\n\
00102 1 0\n\
00103 @end group\n\
00104 @end example\n\
00105 \n\
00106 The matrix is not required to be square.\n\
00107 \n\
00108 When called with two or three output arguments and a spare input matrix,\n\
00109 @code{lu} does not attempt to perform sparsity preserving column\n\
00110 permutations. Called with a fourth output argument, the sparsity\n\
00111 preserving column transformation @var{Q} is returned, such that\n\
00112 @code{@var{P} * @var{A} * @var{Q} = @var{L} * @var{U}}.\n\
00113 \n\
00114 Called with a fifth output argument and a sparse input matrix,\n\
00115 @code{lu} attempts to use a scaling factor @var{R} on the input matrix\n\
00116 such that\n\
00117 @code{@var{P} * (@var{R} \\ @var{A}) * @var{Q} = @var{L} * @var{U}}.\n\
00118 This typically leads to a sparser and more stable factorization.\n\
00119 \n\
00120 An additional input argument @var{thres}, that defines the pivoting\n\
00121 threshold can be given. @var{thres} can be a scalar, in which case\n\
00122 it defines the @sc{umfpack} pivoting tolerance for both symmetric and\n\
00123 unsymmetric cases. If @var{thres} is a 2-element vector, then the first\n\
00124 element defines the pivoting tolerance for the unsymmetric @sc{umfpack}\n\
00125 pivoting strategy and the second for the symmetric strategy. By default,\n\
00126 the values defined by @code{spparms} are used ([0.1, 0.001]).\n\
00127 \n\
00128 Given the string argument 'vector', @code{lu} returns the values of @var{P}\n\
00129 and @var{Q} as vector values, such that for full matrix, @code{@var{A}\n\
00130 (@var{P},:) = @var{L} * @var{U}}, and @code{@var{R}(@var{P},:) * @var{A}\n\
00131 (:, @var{Q}) = @var{L} * @var{U}}.\n\
00132 \n\
00133 With two output arguments, returns the permuted forms of the upper and\n\
00134 lower triangular matrices, such that @code{@var{A} = @var{L} * @var{U}}.\n\
00135 With one output argument @var{y}, then the matrix returned by the @sc{lapack}\n\
00136 routines is returned. If the input matrix is sparse then the matrix @var{L}\n\
00137 is embedded into @var{U} to give a return value similar to the full case.\n\
00138 For both full and sparse matrices, @code{lu} loses the permutation\n\
00139 information.\n\
00140 @end deftypefn")
00141 {
00142 octave_value_list retval;
00143 int nargin = args.length ();
00144 bool issparse = (nargin > 0 && args(0).is_sparse_type ());
00145 bool scale = (nargout == 5);
00146
00147 if (nargin < 1 || (issparse && (nargin > 3 || nargout > 5))
00148 || (!issparse && (nargin > 2 || nargout > 3)))
00149 {
00150 print_usage ();
00151 return retval;
00152 }
00153
00154 bool vecout = false;
00155 Matrix thres;
00156
00157 int n = 1;
00158 while (n < nargin && ! error_state)
00159 {
00160 if (args (n).is_string ())
00161 {
00162 std::string tmp = args(n++).string_value ();
00163
00164 if (! error_state )
00165 {
00166 if (tmp.compare ("vector") == 0)
00167 vecout = true;
00168 else
00169 error ("lu: unrecognized string argument");
00170 }
00171 }
00172 else
00173 {
00174 Matrix tmp = args(n++).matrix_value ();
00175
00176 if (! error_state )
00177 {
00178 if (!issparse)
00179 error ("lu: can not define pivoting threshold THRES for full matrices");
00180 else if (tmp.nelem () == 1)
00181 {
00182 thres.resize(1,2);
00183 thres(0) = tmp(0);
00184 thres(1) = tmp(0);
00185 }
00186 else if (tmp.nelem () == 2)
00187 thres = tmp;
00188 else
00189 error ("lu: expecting 2-element vector for THRES");
00190 }
00191 }
00192 }
00193
00194 octave_value arg = args(0);
00195
00196 octave_idx_type nr = arg.rows ();
00197 octave_idx_type nc = arg.columns ();
00198
00199 int arg_is_empty = empty_arg ("lu", nr, nc);
00200
00201 if (issparse)
00202 {
00203 if (arg_is_empty < 0)
00204 return retval;
00205 else if (arg_is_empty > 0)
00206 return octave_value_list (5, SparseMatrix ());
00207
00208 ColumnVector Qinit;
00209 if (nargout < 4)
00210 {
00211 Qinit.resize (nc);
00212 for (octave_idx_type i = 0; i < nc; i++)
00213 Qinit (i) = i;
00214 }
00215
00216 if (arg.is_real_type ())
00217 {
00218 SparseMatrix m = arg.sparse_matrix_value ();
00219
00220 switch (nargout)
00221 {
00222 case 0:
00223 case 1:
00224 case 2:
00225 {
00226 SparseLU fact (m, Qinit, thres, false, true);
00227
00228 if (nargout < 2)
00229 retval (0) = fact.Y ();
00230 else
00231 {
00232 PermMatrix P = fact.Pr_mat ();
00233 SparseMatrix L = P.transpose () * fact.L ();
00234 retval(1) = octave_value (fact.U (),
00235 MatrixType (MatrixType::Upper));
00236
00237 retval(0) = octave_value (L,
00238 MatrixType (MatrixType::Permuted_Lower,
00239 nr, fact.row_perm ()));
00240 }
00241 }
00242 break;
00243
00244 case 3:
00245 {
00246 SparseLU fact (m, Qinit, thres, false, true);
00247
00248 if (vecout)
00249 retval (2) = fact.Pr_vec ();
00250 else
00251 retval(2) = fact.Pr_mat ();
00252
00253 retval(1) = octave_value (fact.U (),
00254 MatrixType (MatrixType::Upper));
00255 retval(0) = octave_value (fact.L (),
00256 MatrixType (MatrixType::Lower));
00257 }
00258 break;
00259
00260 case 4:
00261 default:
00262 {
00263 SparseLU fact (m, thres, scale);
00264
00265 if (scale)
00266 retval(4) = fact.R ();
00267
00268 if (vecout)
00269 {
00270 retval(3) = fact.Pc_vec ();
00271 retval(2) = fact.Pr_vec ();
00272 }
00273 else
00274 {
00275 retval(3) = fact.Pc_mat ();
00276 retval(2) = fact.Pr_mat ();
00277 }
00278 retval(1) = octave_value (fact.U (),
00279 MatrixType (MatrixType::Upper));
00280 retval(0) = octave_value (fact.L (),
00281 MatrixType (MatrixType::Lower));
00282 }
00283 break;
00284 }
00285 }
00286 else if (arg.is_complex_type ())
00287 {
00288 SparseComplexMatrix m = arg.sparse_complex_matrix_value ();
00289
00290 switch (nargout)
00291 {
00292 case 0:
00293 case 1:
00294 case 2:
00295 {
00296 SparseComplexLU fact (m, Qinit, thres, false, true);
00297
00298 if (nargout < 2)
00299 retval (0) = fact.Y ();
00300 else
00301 {
00302 PermMatrix P = fact.Pr_mat ();
00303 SparseComplexMatrix L = P.transpose () * fact.L ();
00304 retval(1) = octave_value (fact.U (),
00305 MatrixType (MatrixType::Upper));
00306
00307 retval(0) = octave_value (L,
00308 MatrixType (MatrixType::Permuted_Lower,
00309 nr, fact.row_perm ()));
00310 }
00311 }
00312 break;
00313
00314 case 3:
00315 {
00316 SparseComplexLU fact (m, Qinit, thres, false, true);
00317
00318 if (vecout)
00319 retval (2) = fact.Pr_vec ();
00320 else
00321 retval(2) = fact.Pr_mat ();
00322
00323 retval(1) = octave_value (fact.U (),
00324 MatrixType (MatrixType::Upper));
00325 retval(0) = octave_value (fact.L (),
00326 MatrixType (MatrixType::Lower));
00327 }
00328 break;
00329
00330 case 4:
00331 default:
00332 {
00333 SparseComplexLU fact (m, thres, scale);
00334
00335 if (scale)
00336 retval(4) = fact.R ();
00337
00338 if (vecout)
00339 {
00340 retval(3) = fact.Pc_vec ();
00341 retval(2) = fact.Pr_vec ();
00342 }
00343 else
00344 {
00345 retval(3) = fact.Pc_mat ();
00346 retval(2) = fact.Pr_mat ();
00347 }
00348 retval(1) = octave_value (fact.U (),
00349 MatrixType (MatrixType::Upper));
00350 retval(0) = octave_value (fact.L (),
00351 MatrixType (MatrixType::Lower));
00352 }
00353 break;
00354 }
00355 }
00356 else
00357 gripe_wrong_type_arg ("lu", arg);
00358 }
00359 else
00360 {
00361 if (arg_is_empty < 0)
00362 return retval;
00363 else if (arg_is_empty > 0)
00364 return octave_value_list (3, Matrix ());
00365
00366 if (arg.is_real_type ())
00367 {
00368 if (arg.is_single_type ())
00369 {
00370 FloatMatrix m = arg.float_matrix_value ();
00371
00372 if (! error_state)
00373 {
00374 FloatLU fact (m);
00375
00376 switch (nargout)
00377 {
00378 case 0:
00379 case 1:
00380 retval(0) = fact.Y ();
00381 break;
00382
00383 case 2:
00384 {
00385 PermMatrix P = fact.P ();
00386 FloatMatrix L = P.transpose () * fact.L ();
00387 retval(1) = get_lu_u (fact);
00388 retval(0) = L;
00389 }
00390 break;
00391
00392 case 3:
00393 default:
00394 {
00395 if (vecout)
00396 retval(2) = fact.P_vec ();
00397 else
00398 retval(2) = fact.P ();
00399 retval(1) = get_lu_u (fact);
00400 retval(0) = get_lu_l (fact);
00401 }
00402 break;
00403 }
00404 }
00405 }
00406 else
00407 {
00408 Matrix m = arg.matrix_value ();
00409
00410 if (! error_state)
00411 {
00412 LU fact (m);
00413
00414 switch (nargout)
00415 {
00416 case 0:
00417 case 1:
00418 retval(0) = fact.Y ();
00419 break;
00420
00421 case 2:
00422 {
00423 PermMatrix P = fact.P ();
00424 Matrix L = P.transpose () * fact.L ();
00425 retval(1) = get_lu_u (fact);
00426 retval(0) = L;
00427 }
00428 break;
00429
00430 case 3:
00431 default:
00432 {
00433 if (vecout)
00434 retval(2) = fact.P_vec ();
00435 else
00436 retval(2) = fact.P ();
00437 retval(1) = get_lu_u (fact);
00438 retval(0) = get_lu_l (fact);
00439 }
00440 break;
00441 }
00442 }
00443 }
00444 }
00445 else if (arg.is_complex_type ())
00446 {
00447 if (arg.is_single_type ())
00448 {
00449 FloatComplexMatrix m = arg.float_complex_matrix_value ();
00450
00451 if (! error_state)
00452 {
00453 FloatComplexLU fact (m);
00454
00455 switch (nargout)
00456 {
00457 case 0:
00458 case 1:
00459 retval(0) = fact.Y ();
00460 break;
00461
00462 case 2:
00463 {
00464 PermMatrix P = fact.P ();
00465 FloatComplexMatrix L = P.transpose () * fact.L ();
00466 retval(1) = get_lu_u (fact);
00467 retval(0) = L;
00468 }
00469 break;
00470
00471 case 3:
00472 default:
00473 {
00474 if (vecout)
00475 retval(2) = fact.P_vec ();
00476 else
00477 retval(2) = fact.P ();
00478 retval(1) = get_lu_u (fact);
00479 retval(0) = get_lu_l (fact);
00480 }
00481 break;
00482 }
00483 }
00484 }
00485 else
00486 {
00487 ComplexMatrix m = arg.complex_matrix_value ();
00488
00489 if (! error_state)
00490 {
00491 ComplexLU fact (m);
00492
00493 switch (nargout)
00494 {
00495 case 0:
00496 case 1:
00497 retval(0) = fact.Y ();
00498 break;
00499
00500 case 2:
00501 {
00502 PermMatrix P = fact.P ();
00503 ComplexMatrix L = P.transpose () * fact.L ();
00504 retval(1) = get_lu_u (fact);
00505 retval(0) = L;
00506 }
00507 break;
00508
00509 case 3:
00510 default:
00511 {
00512 if (vecout)
00513 retval(2) = fact.P_vec ();
00514 else
00515 retval(2) = fact.P ();
00516 retval(1) = get_lu_u (fact);
00517 retval(0) = get_lu_l (fact);
00518 }
00519 break;
00520 }
00521 }
00522 }
00523 }
00524 else
00525 gripe_wrong_type_arg ("lu", arg);
00526 }
00527
00528 return retval;
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 static
00589 bool check_lu_dims (const octave_value& l, const octave_value& u,
00590 const octave_value& p)
00591 {
00592 octave_idx_type m = l.rows (), k = u.rows (), n = u.columns ();
00593 return ((l.ndims () == 2 && u.ndims () == 2 && k == l.columns ())
00594 && k == std::min (m, n) &&
00595 (p.is_undefined () || p.rows () == m));
00596 }
00597
00598 DEFUN_DLD (luupdate, args, ,
00599 "-*- texinfo -*-\n\
00600 @deftypefn {Loadable Function} {[@var{L}, @var{U}] =} luupdate (@var{L}, @var{U}, @var{x}, @var{y})\n\
00601 @deftypefnx {Loadable Function} {[@var{L}, @var{U}, @var{P}] =} luupdate (@var{L}, @var{U}, @var{P}, @var{x}, @var{y})\n\
00602 Given an LU@tie{}factorization of a real or complex matrix\n\
00603 @w{@var{A} = @var{L}*@var{U}}, @var{L}@tie{}lower unit trapezoidal and\n\
00604 @var{U}@tie{}upper trapezoidal, return the LU@tie{}factorization\n\
00605 of @w{@var{A} + @var{x}*@var{y}.'}, where @var{x} and @var{y} are\n\
00606 column vectors (rank-1 update) or matrices with equal number of columns\n\
00607 (rank-k update).\n\
00608 Optionally, row-pivoted updating can be used by supplying\n\
00609 a row permutation (pivoting) matrix @var{P};\n\
00610 in that case, an updated permutation matrix is returned.\n\
00611 Note that if @var{L}, @var{U}, @var{P} is a pivoted LU@tie{}factorization\n\
00612 as obtained by @code{lu}:\n\
00613 \n\
00614 @example\n\
00615 [@var{L}, @var{U}, @var{P}] = lu (@var{A});\n\
00616 @end example\n\
00617 \n\
00618 @noindent\n\
00619 then a factorization of @xcode{@var{A}+@var{x}*@var{y}.'} can be obtained\n\
00620 either as\n\
00621 \n\
00622 @example\n\
00623 [@var{L1}, @var{U1}] = lu (@var{L}, @var{U}, @var{P}*@var{x}, @var{y})\n\
00624 @end example\n\
00625 \n\
00626 @noindent\n\
00627 or\n\
00628 \n\
00629 @example\n\
00630 [@var{L1}, @var{U1}, @var{P1}] = lu (@var{L}, @var{U}, @var{P}, @var{x}, @var{y})\n\
00631 @end example\n\
00632 \n\
00633 The first form uses the unpivoted algorithm, which is faster, but less\n\
00634 stable. The second form uses a slower pivoted algorithm, which is more\n\
00635 stable.\n\
00636 \n\
00637 The matrix case is done as a sequence of rank-1 updates;\n\
00638 thus, for large enough k, it will be both faster and more accurate to\n\
00639 recompute the factorization from scratch.\n\
00640 @seealso{lu, qrupdate, cholupdate}\n\
00641 @end deftypefn")
00642 {
00643 octave_idx_type nargin = args.length ();
00644 octave_value_list retval;
00645
00646 bool pivoted = nargin == 5;
00647
00648 if (nargin != 4 && nargin != 5)
00649 {
00650 print_usage ();
00651 return retval;
00652 }
00653
00654 octave_value argl = args(0);
00655 octave_value argu = args(1);
00656 octave_value argp = pivoted ? args(2) : octave_value ();
00657 octave_value argx = args(2 + pivoted);
00658 octave_value argy = args(3 + pivoted);
00659
00660 if (argl.is_numeric_type () && argu.is_numeric_type ()
00661 && argx.is_numeric_type () && argy.is_numeric_type ()
00662 && (! pivoted || argp.is_perm_matrix ()))
00663 {
00664 if (check_lu_dims (argl, argu, argp))
00665 {
00666 PermMatrix P = (pivoted
00667 ? argp.perm_matrix_value ()
00668 : PermMatrix::eye (argl.rows ()));
00669
00670 if (argl.is_real_type ()
00671 && argu.is_real_type ()
00672 && argx.is_real_type ()
00673 && argy.is_real_type ())
00674 {
00675
00676 if (argl.is_single_type ()
00677 || argu.is_single_type ()
00678 || argx.is_single_type ()
00679 || argy.is_single_type ())
00680 {
00681 FloatMatrix L = argl.float_matrix_value ();
00682 FloatMatrix U = argu.float_matrix_value ();
00683 FloatMatrix x = argx.float_matrix_value ();
00684 FloatMatrix y = argy.float_matrix_value ();
00685
00686 FloatLU fact (L, U, P);
00687 if (pivoted)
00688 fact.update_piv (x, y);
00689 else
00690 fact.update (x, y);
00691
00692 if (pivoted)
00693 retval(2) = fact.P ();
00694 retval(1) = get_lu_u (fact);
00695 retval(0) = get_lu_l (fact);
00696 }
00697 else
00698 {
00699 Matrix L = argl.matrix_value ();
00700 Matrix U = argu.matrix_value ();
00701 Matrix x = argx.matrix_value ();
00702 Matrix y = argy.matrix_value ();
00703
00704 LU fact (L, U, P);
00705 if (pivoted)
00706 fact.update_piv (x, y);
00707 else
00708 fact.update (x, y);
00709
00710 if (pivoted)
00711 retval(2) = fact.P ();
00712 retval(1) = get_lu_u (fact);
00713 retval(0) = get_lu_l (fact);
00714 }
00715 }
00716 else
00717 {
00718
00719 if (argl.is_single_type ()
00720 || argu.is_single_type ()
00721 || argx.is_single_type ()
00722 || argy.is_single_type ())
00723 {
00724 FloatComplexMatrix L = argl.float_complex_matrix_value ();
00725 FloatComplexMatrix U = argu.float_complex_matrix_value ();
00726 FloatComplexMatrix x = argx.float_complex_matrix_value ();
00727 FloatComplexMatrix y = argy.float_complex_matrix_value ();
00728
00729 FloatComplexLU fact (L, U, P);
00730 if (pivoted)
00731 fact.update_piv (x, y);
00732 else
00733 fact.update (x, y);
00734
00735 if (pivoted)
00736 retval(2) = fact.P ();
00737 retval(1) = get_lu_u (fact);
00738 retval(0) = get_lu_l (fact);
00739 }
00740 else
00741 {
00742 ComplexMatrix L = argl.complex_matrix_value ();
00743 ComplexMatrix U = argu.complex_matrix_value ();
00744 ComplexMatrix x = argx.complex_matrix_value ();
00745 ComplexMatrix y = argy.complex_matrix_value ();
00746
00747 ComplexLU fact (L, U, P);
00748 if (pivoted)
00749 fact.update_piv (x, y);
00750 else
00751 fact.update (x, y);
00752
00753 if (pivoted)
00754 retval(2) = fact.P ();
00755 retval(1) = get_lu_u (fact);
00756 retval(0) = get_lu_l (fact);
00757 }
00758 }
00759 }
00760 else
00761 error ("luupdate: dimension mismatch");
00762 }
00763 else
00764 error ("luupdate: L, U, X, and Y must be numeric");
00765
00766 return retval;
00767 }
00768
00769
00770
00771
00772
00773
00774
00775
00776
00777
00778
00779
00780
00781
00782
00783
00784
00785
00786
00787
00788
00789
00790
00791
00792
00793
00794
00795
00796
00797
00798
00799
00800
00801
00802
00803
00804
00805
00806
00807
00808
00809
00810
00811
00812
00813
00814
00815
00816
00817
00818
00819
00820
00821
00822
00823
00824
00825
00826
00827
00828
00829
00830
00831
00832
00833
00834
00835
00836
00837
00838
00839
00840
00841
00842
00843
00844
00845
00846
00847
00848
00849
00850
00851
00852
00853
00854
00855
00856
00857
00858
00859