Mercurial > octave
view libinterp/corefcn/lu.cc @ 21966:112b20240c87
move docstrings in C++ files out of C strings and into comments
* __contourc__.cc, __dispatch__.cc, __dsearchn__.cc, __ichol__.cc,
__ilu__.cc, __lin_interpn__.cc, __luinc__.cc, __magick_read__.cc,
__pchip_deriv__.cc, __qp__.cc, balance.cc, besselj.cc, betainc.cc,
bitfcns.cc, bsxfun.cc, cellfun.cc, colloc.cc, conv2.cc, daspk.cc,
dasrt.cc, dassl.cc, data.cc, debug.cc, defaults.cc, det.cc, dirfns.cc,
dlmread.cc, dot.cc, eig.cc, ellipj.cc, error.cc, fft.cc, fft2.cc,
fftn.cc, file-io.cc, filter.cc, find.cc, gammainc.cc, gcd.cc,
getgrent.cc, getpwent.cc, getrusage.cc, givens.cc, graphics.cc,
hash.cc, help.cc, hess.cc, hex2num.cc, input.cc, inv.cc, kron.cc,
load-path.cc, load-save.cc, lookup.cc, ls-oct-text.cc, lsode.cc,
lu.cc, mappers.cc, matrix_type.cc, max.cc, mgorth.cc, nproc.cc,
oct-hist.cc, octave-link.cc, ordschur.cc, pager.cc, pinv.cc,
pr-output.cc, profiler.cc, psi.cc, pt-jit.cc, quad.cc, quadcc.cc,
qz.cc, rand.cc, rcond.cc, regexp.cc, schur.cc, sighandlers.cc,
sparse.cc, spparms.cc, sqrtm.cc, str2double.cc, strfind.cc, strfns.cc,
sub2ind.cc, svd.cc, sylvester.cc, symtab.cc, syscalls.cc, sysdep.cc,
time.cc, toplev.cc, tril.cc, tsearch.cc, typecast.cc, urlwrite.cc,
utils.cc, variables.cc, __delaunayn__.cc, __eigs__.cc,
__fltk_uigetfile__.cc, __glpk__.cc, __init_fltk__.cc,
__init_gnuplot__.cc, __osmesa_print__.cc, __voronoi__.cc, amd.cc,
audiodevinfo.cc, audioread.cc, ccolamd.cc, chol.cc, colamd.cc,
convhulln.cc, dmperm.cc, fftw.cc, qr.cc, symbfact.cc, symrcm.cc,
ov-base.cc, ov-bool-mat.cc, ov-cell.cc, ov-class.cc, ov-classdef.cc,
ov-fcn-handle.cc, ov-fcn-inline.cc, ov-flt-re-mat.cc, ov-int16.cc,
ov-int32.cc, ov-int64.cc, ov-int8.cc, ov-java.cc, ov-null-mat.cc,
ov-oncleanup.cc, ov-range.cc, ov-re-mat.cc, ov-struct.cc,
ov-typeinfo.cc, ov-uint16.cc, ov-uint32.cc, ov-uint64.cc, ov-uint8.cc,
ov-usr-fcn.cc, ov.cc, octave.cc, pt-arg-list.cc, pt-binop.cc,
pt-eval.cc, pt-mat.cc, lex.ll, oct-parse.in.yy:
Docstrings are now comments instead of C strings.
* build-aux/mk-opts.pl: Emit docstrings as comments instead of C
strings.
* DASPK-opts.in, LSODE-opts.in: Don't quote " in docstring fragments.
* builtins.h: Include builtin-defun-decls.h unconditionally.
* defun.h (DEFUN, DEFUNX, DEFCONSTFUN): Simply emit declaration.
(DEFALIAS): Always expand to nothing.
* defun-dld.h: No special macro expansions for MAKE_BUILTINS.
(DEFUN_DLD): Use FORWARD_DECLARE_FUN.
(DEFUNX_DLD): Use FORWARD_DECLARE_FUNX.
* defun-int.h: No special macro expansions for MAKE_BUILTINS.
(FORWARD_DECLARE_FUN, FORWARD_DECLARE_FUNX): New macros.
(DEFINE_FUN_INSTALLER_FUN): If compiling an Octave source file, pass
"external-doc" to DEFINE_FUNX_INSTALLER_FUN.
(DEFUN_INTERNAL, DEFCONSTFUN_INTERNAL, DEFUNX_INTERNAL,
DEFALIAS_INTERNAL): Delete.
* common.mk (move_if_change_rule): New macro.
(simple_move_if_change_rule): Define using move_if_change_rule.
* find-defun-files.sh (DEFUN_PATTERN): Update. Don't transform file
name extension to ".df".
* libinterp/mk-pkg-add, gendoc.pl: Operate directly on source files.
* mkbuiltins: New argument, SRCDIR. Operate directly on source files.
* mkdefs: Delete.
* libinterp/module.mk (BUILT_SOURCES): Update list to contain only
files included in other source files.
(GENERATED_MAKE_BUILTINS_INCS, DEF_FILES): Delete.
(LIBINTERP_BUILT_DISTFILES): Include $(OPT_HANDLERS) here.
(LIBINTERP_BUILT_NODISTFILES): Not here. Remove $(ALL_DEF_FILES from
the list.
(libinterp_EXTRA_DIST): Remove mkdefs from the list.
(FOUND_DEFUN_FILES): Rename from SRC_DEF_FILES.
(DLDFCN_DEFUN_FILES): Rename from DLDFCN_DEF_FILES.
(SRC_DEFUN_FILES): Rename from SRC_DEF_FILES.
(ALL_DEFUN_FILES): Rename from ALL_DEF_FILES.
(%.df: %.cc): Delete pattern rule.
(libinterp/build-env-features.cc, libinterp/builtins.cc,
libinterp/dldfcn/PKG_ADD): Use mv instead of move-if-change.
(libinterp/builtins.cc, libinterp/builtin-defun-decls.h):
Update mkbuiltins command.
($(srcdir)/libinterp/DOCSTRINGS): Update gendoc.pl command.
* liboctave/module.mk (BUILT_SOURCES): Don't include
liboctave-build-info.cc in the list.
author | John W. Eaton <jwe@octave.org> |
---|---|
date | Tue, 21 Jun 2016 16:07:51 -0400 |
parents | aba2e6293dd8 |
children | 6ca3acf5fad8 |
line wrap: on
line source
/* Copyright (C) 1996-2015 John W. Eaton This file is part of Octave. Octave is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version. Octave is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with Octave; see the file COPYING. If not, see <http://www.gnu.org/licenses/>. */ #if defined (HAVE_CONFIG_H) # include "config.h" #endif #include "lu.h" #include "sparse-lu.h" #include "defun.h" #include "error.h" #include "errwarn.h" #include "ovl.h" #include "utils.h" #include "ov-re-sparse.h" #include "ov-cx-sparse.h" template <typename MT> static octave_value get_lu_l (const lu<MT>& fact) { MT L = fact.L (); if (L.is_square ()) return octave_value (L, MatrixType (MatrixType::Lower)); else return L; } template <typename MT> static octave_value get_lu_u (const lu<MT>& fact) { MT U = fact.U (); if (U.is_square () && fact.regular ()) return octave_value (U, MatrixType (MatrixType::Upper)); else return U; } DEFUN (lu, args, nargout, doc: /* -*- texinfo -*- @deftypefn {} {[@var{L}, @var{U}] =} lu (@var{A}) @deftypefnx {} {[@var{L}, @var{U}, @var{P}] =} lu (@var{A}) @deftypefnx {} {[@var{L}, @var{U}, @var{P}, @var{Q}] =} lu (@var{S}) @deftypefnx {} {[@var{L}, @var{U}, @var{P}, @var{Q}, @var{R}] =} lu (@var{S}) @deftypefnx {} {[@dots{}] =} lu (@var{S}, @var{thres}) @deftypefnx {} {@var{y} =} lu (@dots{}) @deftypefnx {} {[@dots{}] =} lu (@dots{}, "vector") @cindex LU decomposition Compute the LU@tie{}decomposition of @var{A}. If @var{A} is full then subroutines from @sc{lapack} are used, and if @var{A} is sparse then @sc{umfpack} is used. The result is returned in a permuted form, according to the optional return value @var{P}. For example, given the matrix @code{a = [1, 2; 3, 4]}, @example [l, u, p] = lu (@var{a}) @end example @noindent returns @example @group l = 1.00000 0.00000 0.33333 1.00000 u = 3.00000 4.00000 0.00000 0.66667 p = 0 1 1 0 @end group @end example The matrix is not required to be square. When called with two or three output arguments and a sparse input matrix, @code{lu} does not attempt to perform sparsity preserving column permutations. Called with a fourth output argument, the sparsity preserving column transformation @var{Q} is returned, such that @code{@var{P} * @var{A} * @var{Q} = @var{L} * @var{U}}. Called with a fifth output argument and a sparse input matrix, @code{lu} attempts to use a scaling factor @var{R} on the input matrix such that @code{@var{P} * (@var{R} \ @var{A}) * @var{Q} = @var{L} * @var{U}}. This typically leads to a sparser and more stable factorization. An additional input argument @var{thres}, that defines the pivoting threshold can be given. @var{thres} can be a scalar, in which case it defines the @sc{umfpack} pivoting tolerance for both symmetric and unsymmetric cases. If @var{thres} is a 2-element vector, then the first element defines the pivoting tolerance for the unsymmetric @sc{umfpack} pivoting strategy and the second for the symmetric strategy. By default, the values defined by @code{spparms} are used ([0.1, 0.001]). Given the string argument @qcode{"vector"}, @code{lu} returns the values of @var{P} and @var{Q} as vector values, such that for full matrix, @code{@var{A}(@var{P},:) = @var{L} * @var{U}}, and @code{@var{R}(@var{P},:) * @var{A}(:,@var{Q}) = @var{L} * @var{U}}. With two output arguments, returns the permuted forms of the upper and lower triangular matrices, such that @code{@var{A} = @var{L} * @var{U}}. With one output argument @var{y}, then the matrix returned by the @sc{lapack} routines is returned. If the input matrix is sparse then the matrix @var{L} is embedded into @var{U} to give a return value similar to the full case. For both full and sparse matrices, @code{lu} loses the permutation information. @seealso{luupdate, ilu, chol, hess, qr, qz, schur, svd} @end deftypefn */) { int nargin = args.length (); bool issparse = (nargin > 0 && args(0).is_sparse_type ()); if (nargin < 1 || (issparse && nargin > 3) || (! issparse && nargin > 2)) print_usage (); bool vecout = false; Matrix thres; int n = 1; while (n < nargin) { if (args(n).is_string ()) { std::string tmp = args(n++).string_value (); if (tmp == "vector") vecout = true; else error ("lu: unrecognized string argument"); } else { if (! issparse) error ("lu: can not define pivoting threshold THRES for full matrices"); Matrix tmp = args(n++).matrix_value (); if (tmp.numel () == 1) { thres.resize (1,2); thres(0) = tmp(0); thres(1) = tmp(0); } else if (tmp.numel () == 2) thres = tmp; else error ("lu: THRES must be a 1 or 2-element vector"); } } octave_value_list retval; bool scale = (nargout == 5); octave_value arg = args(0); octave_idx_type nr = arg.rows (); octave_idx_type nc = arg.columns (); int arg_is_empty = empty_arg ("lu", nr, nc); if (issparse) { if (arg_is_empty < 0) return ovl (); else if (arg_is_empty > 0) return octave_value_list (5, SparseMatrix ()); if (arg.is_real_type ()) { SparseMatrix m = arg.sparse_matrix_value (); if (nargout < 4) { ColumnVector Qinit (nc); for (octave_idx_type i = 0; i < nc; i++) Qinit(i) = i; sparse_lu<SparseMatrix> fact (m, Qinit, thres, false, true); if (nargout < 2) retval(0) = fact.Y (); else { retval.resize (nargout == 3 ? 3 : 2); retval(1) = octave_value ( fact.U () * fact.Pc_mat ().transpose (), MatrixType (MatrixType::Permuted_Upper, nc, fact.col_perm ())); PermMatrix P = fact.Pr_mat (); SparseMatrix L = fact.L (); if (nargout == 2) retval(0) = octave_value (P.transpose () * L, MatrixType (MatrixType::Permuted_Lower, nr, fact.row_perm ())); else { retval(0) = L; if (vecout) retval(2) = fact.Pr_vec(); else retval(2) = P; } } } else { retval.resize (scale ? 5 : 4); sparse_lu<SparseMatrix> fact (m, thres, scale); retval(0) = octave_value (fact.L (), MatrixType (MatrixType::Lower)); retval(1) = octave_value (fact.U (), MatrixType (MatrixType::Upper)); if (vecout) { retval(2) = fact.Pr_vec (); retval(3) = fact.Pc_vec (); } else { retval(2) = fact.Pr_mat (); retval(3) = fact.Pc_mat (); } if (scale) retval(4) = fact.R (); } } else if (arg.is_complex_type ()) { SparseComplexMatrix m = arg.sparse_complex_matrix_value (); if (nargout < 4) { ColumnVector Qinit (nc); for (octave_idx_type i = 0; i < nc; i++) Qinit(i) = i; sparse_lu<SparseComplexMatrix> fact (m, Qinit, thres, false, true); if (nargout < 2) retval(0) = fact.Y (); else { retval.resize (nargout == 3 ? 3 : 2); retval(1) = octave_value ( fact.U () * fact.Pc_mat ().transpose (), MatrixType (MatrixType::Permuted_Upper, nc, fact.col_perm ())); PermMatrix P = fact.Pr_mat (); SparseComplexMatrix L = fact.L (); if (nargout == 2) retval(0) = octave_value (P.transpose () * L, MatrixType (MatrixType::Permuted_Lower, nr, fact.row_perm ())); else { retval(0) = L; if (vecout) retval(2) = fact.Pr_vec(); else retval(2) = P; } } } else { retval.resize (scale ? 5 : 4); sparse_lu<SparseComplexMatrix> fact (m, thres, scale); retval(0) = octave_value (fact.L (), MatrixType (MatrixType::Lower)); retval(1) = octave_value (fact.U (), MatrixType (MatrixType::Upper)); if (vecout) { retval(2) = fact.Pr_vec (); retval(3) = fact.Pc_vec (); } else { retval(2) = fact.Pr_mat (); retval(3) = fact.Pc_mat (); } if (scale) retval(4) = fact.R (); } } else err_wrong_type_arg ("lu", arg); } else { if (arg_is_empty < 0) return ovl (); else if (arg_is_empty > 0) return octave_value_list (3, Matrix ()); if (arg.is_real_type ()) { if (arg.is_single_type ()) { FloatMatrix m = arg.float_matrix_value (); lu<FloatMatrix> fact (m); switch (nargout) { case 0: case 1: retval = ovl (fact.Y ()); break; case 2: { PermMatrix P = fact.P (); FloatMatrix L = P.transpose () * fact.L (); retval = ovl (L, get_lu_u (fact)); } break; case 3: default: { if (vecout) retval = ovl (get_lu_l (fact), get_lu_u (fact), fact.P_vec ()); else retval = ovl (get_lu_l (fact), get_lu_u (fact), fact.P ()); } break; } } else { Matrix m = arg.matrix_value (); lu<Matrix> fact (m); switch (nargout) { case 0: case 1: retval = ovl (fact.Y ()); break; case 2: { PermMatrix P = fact.P (); Matrix L = P.transpose () * fact.L (); retval = ovl (L, get_lu_u (fact)); } break; case 3: default: { if (vecout) retval = ovl (get_lu_l (fact), get_lu_u (fact), fact.P_vec ()); else retval = ovl (get_lu_l (fact), get_lu_u (fact), fact.P ()); } break; } } } else if (arg.is_complex_type ()) { if (arg.is_single_type ()) { FloatComplexMatrix m = arg.float_complex_matrix_value (); lu<FloatComplexMatrix> fact (m); switch (nargout) { case 0: case 1: retval = ovl (fact.Y ()); break; case 2: { PermMatrix P = fact.P (); FloatComplexMatrix L = P.transpose () * fact.L (); retval = ovl (L, get_lu_u (fact)); } break; case 3: default: { if (vecout) retval = ovl (get_lu_l (fact), get_lu_u (fact), fact.P_vec ()); else retval = ovl (get_lu_l (fact), get_lu_u (fact), fact.P ()); } break; } } else { ComplexMatrix m = arg.complex_matrix_value (); lu<ComplexMatrix> fact (m); switch (nargout) { case 0: case 1: retval = ovl (fact.Y ()); break; case 2: { PermMatrix P = fact.P (); ComplexMatrix L = P.transpose () * fact.L (); retval = ovl (L, get_lu_u (fact)); } break; case 3: default: { if (vecout) retval = ovl (get_lu_l (fact), get_lu_u (fact), fact.P_vec ()); else retval = ovl (get_lu_l (fact), get_lu_u (fact), fact.P ()); } break; } } } else err_wrong_type_arg ("lu", arg); } return retval; } /* %!assert(lu ([1, 2; 3, 4]), [3, 4; 1/3, 2/3], eps) %!test %! [l, u] = lu ([1, 2; 3, 4]); %! assert (l, [1/3, 1; 1, 0], sqrt (eps)); %! assert (u, [3, 4; 0, 2/3], sqrt (eps)); %!test %! [l, u, p] = lu ([1, 2; 3, 4]); %! assert (l, [1, 0; 1/3, 1], sqrt (eps)); %! assert (u, [3, 4; 0, 2/3], sqrt (eps)); %! assert (p(:,:), [0, 1; 1, 0], sqrt (eps)); %!test %! [l, u, p] = lu ([1, 2; 3, 4], "vector"); %! assert (l, [1, 0; 1/3, 1], sqrt (eps)); %! assert (u, [3, 4; 0, 2/3], sqrt (eps)); %! assert (p, [2;1], sqrt (eps)); %!test %! [l, u, p] = lu ([1, 2; 3, 4; 5, 6]); %! assert (l, [1, 0; 1/5, 1; 3/5, 1/2], sqrt (eps)); %! assert (u, [5, 6; 0, 4/5], sqrt (eps)); %! assert (p(:,:), [0, 0, 1; 1, 0, 0; 0 1 0], sqrt (eps)); %!assert (lu (single ([1, 2; 3, 4])), single ([3, 4; 1/3, 2/3]), eps ("single")) %!test %! [l, u] = lu (single ([1, 2; 3, 4])); %! assert (l, single ([1/3, 1; 1, 0]), sqrt (eps ("single"))); %! assert (u, single ([3, 4; 0, 2/3]), sqrt (eps ("single"))); %!test %! [l, u, p] = lu (single ([1, 2; 3, 4])); %! assert (l, single ([1, 0; 1/3, 1]), sqrt (eps ("single"))); %! assert (u, single ([3, 4; 0, 2/3]), sqrt (eps ("single"))); %! assert (p(:,:), single ([0, 1; 1, 0]), sqrt (eps ("single"))); %!test %! [l, u, p] = lu (single ([1, 2; 3, 4]), "vector"); %! assert (l, single ([1, 0; 1/3, 1]), sqrt (eps ("single"))); %! assert (u, single ([3, 4; 0, 2/3]), sqrt (eps ("single"))); %! assert (p, single ([2;1]), sqrt (eps ("single"))); %!test %! [l u p] = lu (single ([1, 2; 3, 4; 5, 6])); %! assert (l, single ([1, 0; 1/5, 1; 3/5, 1/2]), sqrt (eps ("single"))); %! assert (u, single ([5, 6; 0, 4/5]), sqrt (eps ("single"))); %! assert (p(:,:), single ([0, 0, 1; 1, 0, 0; 0 1 0]), sqrt (eps ("single"))); %!error lu () %!error <can not define pivoting threshold> lu ([1, 2; 3, 4], 2) %!testif HAVE_UMFPACK %! Bi = [1 2 3 4 5 2 3 6 7 8 4 5 7 8 9]; %! Bj = [1 3 4 5 6 7 8 9 11 12 13 14 15 16 17]; %! Bv = [1 1 1 1 1 1 -1 1 1 1 1 -1 1 -1 1]; %! B = sparse (Bi, Bj, Bv); %! [L, U] = lu (B); %! assert (L*U, B); %! [L, U, P] = lu(B); %! assert (P'*L*U, B); %! [L, U, P, Q] = lu (B); %! assert (P'*L*U*Q', B); */ static bool check_lu_dims (const octave_value& l, const octave_value& u, const octave_value& p) { octave_idx_type m = l.rows (); octave_idx_type k = u.rows (); octave_idx_type n = u.columns (); return ((l.ndims () == 2 && u.ndims () == 2 && k == l.columns ()) && k == std::min (m, n) && (p.is_undefined () || p.rows () == m)); } DEFUN (luupdate, args, , doc: /* -*- texinfo -*- @deftypefn {} {[@var{L}, @var{U}] =} luupdate (@var{L}, @var{U}, @var{x}, @var{y}) @deftypefnx {} {[@var{L}, @var{U}, @var{P}] =} luupdate (@var{L}, @var{U}, @var{P}, @var{x}, @var{y}) Given an LU@tie{}factorization of a real or complex matrix @w{@var{A} = @var{L}*@var{U}}, @var{L}@tie{}lower unit trapezoidal and @var{U}@tie{}upper trapezoidal, return the LU@tie{}factorization of @w{@var{A} + @var{x}*@var{y}.'}, where @var{x} and @var{y} are column vectors (rank-1 update) or matrices with equal number of columns (rank-k update). Optionally, row-pivoted updating can be used by supplying a row permutation (pivoting) matrix @var{P}; in that case, an updated permutation matrix is returned. Note that if @var{L}, @var{U}, @var{P} is a pivoted LU@tie{}factorization as obtained by @code{lu}: @example [@var{L}, @var{U}, @var{P}] = lu (@var{A}); @end example @noindent then a factorization of @tcode{@var{A}+@var{x}*@var{y}.'} can be obtained either as @example [@var{L1}, @var{U1}] = lu (@var{L}, @var{U}, @var{P}*@var{x}, @var{y}) @end example @noindent or @example [@var{L1}, @var{U1}, @var{P1}] = lu (@var{L}, @var{U}, @var{P}, @var{x}, @var{y}) @end example The first form uses the unpivoted algorithm, which is faster, but less stable. The second form uses a slower pivoted algorithm, which is more stable. The matrix case is done as a sequence of rank-1 updates; thus, for large enough k, it will be both faster and more accurate to recompute the factorization from scratch. @seealso{lu, cholupdate, qrupdate} @end deftypefn */) { int nargin = args.length (); if (nargin != 4 && nargin != 5) print_usage (); bool pivoted = (nargin == 5); octave_value argl = args(0); octave_value argu = args(1); octave_value argp = pivoted ? args(2) : octave_value (); octave_value argx = args(2 + pivoted); octave_value argy = args(3 + pivoted); if (! (argl.is_numeric_type () && argu.is_numeric_type () && argx.is_numeric_type () && argy.is_numeric_type () && (! pivoted || argp.is_perm_matrix ()))) error ("luupdate: L, U, X, and Y must be numeric"); if (! check_lu_dims (argl, argu, argp)) error ("luupdate: dimension mismatch"); PermMatrix P = (pivoted ? argp.perm_matrix_value () : PermMatrix::eye (argl.rows ())); if (argl.is_real_type () && argu.is_real_type () && argx.is_real_type () && argy.is_real_type ()) { // all real case if (argl.is_single_type () || argu.is_single_type () || argx.is_single_type () || argy.is_single_type ()) { FloatMatrix L = argl.float_matrix_value (); FloatMatrix U = argu.float_matrix_value (); FloatMatrix x = argx.float_matrix_value (); FloatMatrix y = argy.float_matrix_value (); lu<FloatMatrix> fact (L, U, P); if (pivoted) fact.update_piv (x, y); else fact.update (x, y); if (pivoted) return ovl (get_lu_l (fact), get_lu_u (fact), fact.P ()); else return ovl (get_lu_l (fact), get_lu_u (fact)); } else { Matrix L = argl.matrix_value (); Matrix U = argu.matrix_value (); Matrix x = argx.matrix_value (); Matrix y = argy.matrix_value (); lu<Matrix> fact (L, U, P); if (pivoted) fact.update_piv (x, y); else fact.update (x, y); if (pivoted) return ovl (get_lu_l (fact), get_lu_u (fact), fact.P ()); else return ovl (get_lu_l (fact), get_lu_u (fact)); } } else { // complex case if (argl.is_single_type () || argu.is_single_type () || argx.is_single_type () || argy.is_single_type ()) { FloatComplexMatrix L = argl.float_complex_matrix_value (); FloatComplexMatrix U = argu.float_complex_matrix_value (); FloatComplexMatrix x = argx.float_complex_matrix_value (); FloatComplexMatrix y = argy.float_complex_matrix_value (); lu<FloatComplexMatrix> fact (L, U, P); if (pivoted) fact.update_piv (x, y); else fact.update (x, y); if (pivoted) return ovl (get_lu_l (fact), get_lu_u (fact), fact.P ()); else return ovl (get_lu_l (fact), get_lu_u (fact)); } else { ComplexMatrix L = argl.complex_matrix_value (); ComplexMatrix U = argu.complex_matrix_value (); ComplexMatrix x = argx.complex_matrix_value (); ComplexMatrix y = argy.complex_matrix_value (); lu<ComplexMatrix> fact (L, U, P); if (pivoted) fact.update_piv (x, y); else fact.update (x, y); if (pivoted) return ovl (get_lu_l (fact), get_lu_u (fact), fact.P ()); else return ovl (get_lu_l (fact), get_lu_u (fact)); } } } /* %!shared A, u, v, Ac, uc, vc %! A = [0.091364 0.613038 0.999083; %! 0.594638 0.425302 0.603537; %! 0.383594 0.291238 0.085574; %! 0.265712 0.268003 0.238409; %! 0.669966 0.743851 0.445057 ]; %! %! u = [0.85082; %! 0.76426; %! 0.42883; %! 0.53010; %! 0.80683 ]; %! %! v = [0.98810; %! 0.24295; %! 0.43167 ]; %! %! Ac = [0.620405 + 0.956953i 0.480013 + 0.048806i 0.402627 + 0.338171i; %! 0.589077 + 0.658457i 0.013205 + 0.279323i 0.229284 + 0.721929i; %! 0.092758 + 0.345687i 0.928679 + 0.241052i 0.764536 + 0.832406i; %! 0.912098 + 0.721024i 0.049018 + 0.269452i 0.730029 + 0.796517i; %! 0.112849 + 0.603871i 0.486352 + 0.142337i 0.355646 + 0.151496i ]; %! %! uc = [0.20351 + 0.05401i; %! 0.13141 + 0.43708i; %! 0.29808 + 0.08789i; %! 0.69821 + 0.38844i; %! 0.74871 + 0.25821i ]; %! %! vc = [0.85839 + 0.29468i; %! 0.20820 + 0.93090i; %! 0.86184 + 0.34689i ]; %! %!testif HAVE_QRUPDATE_LUU %! [L,U,P] = lu (A); %! [L,U] = luupdate (L,U,P*u,v); %! assert (norm (vec (tril (L)-L), Inf) == 0); %! assert (norm (vec (triu (U)-U), Inf) == 0); %! assert (norm (vec (P'*L*U - A - u*v.'), Inf) < norm (A)*1e1*eps); %! %!testif HAVE_QRUPDATE_LUU %! [L,U,P] = lu (Ac); %! [L,U] = luupdate (L,U,P*uc,vc); %! assert (norm (vec (tril (L)-L), Inf) == 0); %! assert (norm (vec (triu (U)-U), Inf) == 0); %! assert (norm (vec (P'*L*U - Ac - uc*vc.'), Inf) < norm (Ac)*1e1*eps); %!testif HAVE_QRUPDATE_LUU %! [L,U,P] = lu (single (A)); %! [L,U] = luupdate (L,U,P*single (u), single (v)); %! assert (norm (vec (tril (L)-L), Inf) == 0); %! assert (norm (vec (triu (U)-U), Inf) == 0); %! assert (norm (vec (P'*L*U - single (A) - single (u)*single (v).'), Inf) < norm (single (A))*1e1*eps ("single")); %! %!testif HAVE_QRUPDATE_LUU %! [L,U,P] = lu (single (Ac)); %! [L,U] = luupdate (L,U,P*single (uc),single (vc)); %! assert (norm (vec (tril (L)-L), Inf) == 0); %! assert (norm (vec (triu (U)-U), Inf) == 0); %! assert (norm (vec (P'*L*U - single (Ac) - single (uc)*single (vc).'), Inf) < norm (single (Ac))*1e1*eps ("single")); %!testif HAVE_QRUPDATE_LUU %! [L,U,P] = lu (A); %! [L,U,P] = luupdate (L,U,P,u,v); %! assert (norm (vec (tril (L)-L), Inf) == 0); %! assert (norm (vec (triu (U)-U), Inf) == 0); %! assert (norm (vec (P'*L*U - A - u*v.'), Inf) < norm (A)*1e1*eps); %! %!testif HAVE_QRUPDATE_LUU %! [L,U,P] = lu (A); %! [~,ordcols] = max (P,[],1); %! [~,ordrows] = max (P,[],2); %! P1 = eye (size(P))(:,ordcols); %! P2 = eye (size(P))(ordrows,:); %! assert(P1 == P); %! assert(P2 == P); %! [L,U,P] = luupdate (L,U,P,u,v); %! [L,U,P1] = luupdate (L,U,P1,u,v); %! [L,U,P2] = luupdate (L,U,P2,u,v); %! assert(P1 == P); %! assert(P2 == P); %! %!testif HAVE_QRUPDATE_LUU %! [L,U,P] = lu (Ac); %! [L,U,P] = luupdate (L,U,P,uc,vc); %! assert (norm (vec (tril (L)-L), Inf) == 0); %! assert (norm (vec (triu (U)-U), Inf) == 0); %! assert (norm (vec (P'*L*U - Ac - uc*vc.'), Inf) < norm (Ac)*1e1*eps); %!testif HAVE_QRUPDATE_LUU %! [L,U,P] = lu (single (A)); %! [L,U,P] = luupdate (L,U,P,single (u),single (v)); %! assert (norm (vec (tril (L)-L), Inf) == 0); %! assert (norm (vec (triu (U)-U), Inf) == 0); %! assert (norm (vec (P'*L*U - single (A) - single (u)*single (v).'), Inf) < norm (single (A))*1e1*eps ("single")); %! %!testif HAVE_QRUPDATE_LUU %! [L,U,P] = lu (single (Ac)); %! [L,U,P] = luupdate (L,U,P,single (uc),single (vc)); %! assert (norm (vec (tril (L)-L), Inf) == 0); %! assert (norm (vec (triu (U)-U), Inf) == 0); %! assert (norm (vec (P'*L*U - single (Ac) - single (uc)*single (vc).'), Inf) < norm (single (Ac))*1e1*eps ("single")); */