Mercurial > octave
annotate 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 |
rev | line source |
---|---|
2928 | 1 /* |
2 | |
19697
4197fc428c7d
maint: Update copyright notices for 2015.
John W. Eaton <jwe@octave.org>
parents:
19596
diff
changeset
|
3 Copyright (C) 1996-2015 John W. Eaton |
2928 | 4 |
5 This file is part of Octave. | |
6 | |
7 Octave is free software; you can redistribute it and/or modify it | |
8 under the terms of the GNU General Public License as published by the | |
7016 | 9 Free Software Foundation; either version 3 of the License, or (at your |
10 option) any later version. | |
2928 | 11 |
12 Octave is distributed in the hope that it will be useful, but WITHOUT | |
13 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or | |
14 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License | |
15 for more details. | |
16 | |
17 You should have received a copy of the GNU General Public License | |
7016 | 18 along with Octave; see the file COPYING. If not, see |
19 <http://www.gnu.org/licenses/>. | |
2928 | 20 |
21 */ | |
22 | |
21724
aba2e6293dd8
use "#if ..." consistently instead of "#ifdef" and "#ifndef"
John W. Eaton <jwe@octave.org>
parents:
21547
diff
changeset
|
23 #if defined (HAVE_CONFIG_H) |
21301
40de9f8f23a6
Use '#include "config.h"' rather than <config.h>.
Rik <rik@octave.org>
parents:
21271
diff
changeset
|
24 # include "config.h" |
2928 | 25 #endif |
26 | |
21271
7e67c7f82fc1
better use of templates for lu factorization classes
John W. Eaton <jwe@octave.org>
parents:
21200
diff
changeset
|
27 #include "lu.h" |
21146
ea9c05014809
revamp sparse LU factorization classes
John W. Eaton <jwe@octave.org>
parents:
21139
diff
changeset
|
28 #include "sparse-lu.h" |
2928 | 29 |
15039
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14854
diff
changeset
|
30 #include "defun.h" |
2928 | 31 #include "error.h" |
21100
e39e05d90788
Switch gripe_XXX to either err_XXX or warn_XXX naming scheme.
Rik <rik@octave.org>
parents:
21055
diff
changeset
|
32 #include "errwarn.h" |
20940
48b2ad5ee801
maint: Rename oct-obj.[cc|h] to ovl.[cc|h] for clarity.
Rik <rik@octave.org>
parents:
20921
diff
changeset
|
33 #include "ovl.h" |
2928 | 34 #include "utils.h" |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
35 #include "ov-re-sparse.h" |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
36 #include "ov-cx-sparse.h" |
2928 | 37 |
21139
538b57866b90
consistently use "typename" intead of "class" in template declarations
John W. Eaton <jwe@octave.org>
parents:
21100
diff
changeset
|
38 template <typename MT> |
8869
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
39 static octave_value |
21271
7e67c7f82fc1
better use of templates for lu factorization classes
John W. Eaton <jwe@octave.org>
parents:
21200
diff
changeset
|
40 get_lu_l (const lu<MT>& fact) |
8869
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
41 { |
9715
9f27172fbd1e
auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents:
9709
diff
changeset
|
42 MT L = fact.L (); |
9f27172fbd1e
auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents:
9709
diff
changeset
|
43 if (L.is_square ()) |
9f27172fbd1e
auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents:
9709
diff
changeset
|
44 return octave_value (L, MatrixType (MatrixType::Lower)); |
8869
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
45 else |
9715
9f27172fbd1e
auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents:
9709
diff
changeset
|
46 return L; |
9f27172fbd1e
auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents:
9709
diff
changeset
|
47 } |
8869
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
48 |
21139
538b57866b90
consistently use "typename" intead of "class" in template declarations
John W. Eaton <jwe@octave.org>
parents:
21100
diff
changeset
|
49 template <typename MT> |
9715
9f27172fbd1e
auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents:
9709
diff
changeset
|
50 static octave_value |
21271
7e67c7f82fc1
better use of templates for lu factorization classes
John W. Eaton <jwe@octave.org>
parents:
21200
diff
changeset
|
51 get_lu_u (const lu<MT>& fact) |
9715
9f27172fbd1e
auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents:
9709
diff
changeset
|
52 { |
9f27172fbd1e
auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents:
9709
diff
changeset
|
53 MT U = fact.U (); |
9f27172fbd1e
auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents:
9709
diff
changeset
|
54 if (U.is_square () && fact.regular ()) |
9f27172fbd1e
auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents:
9709
diff
changeset
|
55 return octave_value (U, MatrixType (MatrixType::Upper)); |
9f27172fbd1e
auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents:
9709
diff
changeset
|
56 else |
9f27172fbd1e
auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents:
9709
diff
changeset
|
57 return U; |
8869
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
58 } |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
59 |
15039
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14854
diff
changeset
|
60 DEFUN (lu, args, nargout, |
21966
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
61 doc: /* -*- texinfo -*- |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
62 @deftypefn {} {[@var{L}, @var{U}] =} lu (@var{A}) |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
63 @deftypefnx {} {[@var{L}, @var{U}, @var{P}] =} lu (@var{A}) |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
64 @deftypefnx {} {[@var{L}, @var{U}, @var{P}, @var{Q}] =} lu (@var{S}) |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
65 @deftypefnx {} {[@var{L}, @var{U}, @var{P}, @var{Q}, @var{R}] =} lu (@var{S}) |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
66 @deftypefnx {} {[@dots{}] =} lu (@var{S}, @var{thres}) |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
67 @deftypefnx {} {@var{y} =} lu (@dots{}) |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
68 @deftypefnx {} {[@dots{}] =} lu (@dots{}, "vector") |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
69 @cindex LU decomposition |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
70 Compute the LU@tie{}decomposition of @var{A}. |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
71 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
72 If @var{A} is full then subroutines from @sc{lapack} are used, and if |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
73 @var{A} is sparse then @sc{umfpack} is used. |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
74 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
75 The result is returned in a permuted form, according to the optional return |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
76 value @var{P}. For example, given the matrix @code{a = [1, 2; 3, 4]}, |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
77 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
78 @example |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
79 [l, u, p] = lu (@var{a}) |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
80 @end example |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
81 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
82 @noindent |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
83 returns |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
84 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
85 @example |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
86 @group |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
87 l = |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
88 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
89 1.00000 0.00000 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
90 0.33333 1.00000 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
91 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
92 u = |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
93 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
94 3.00000 4.00000 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
95 0.00000 0.66667 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
96 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
97 p = |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
98 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
99 0 1 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
100 1 0 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
101 @end group |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
102 @end example |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
103 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
104 The matrix is not required to be square. |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
105 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
106 When called with two or three output arguments and a sparse input matrix, |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
107 @code{lu} does not attempt to perform sparsity preserving column |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
108 permutations. Called with a fourth output argument, the sparsity |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
109 preserving column transformation @var{Q} is returned, such that |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
110 @code{@var{P} * @var{A} * @var{Q} = @var{L} * @var{U}}. |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
111 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
112 Called with a fifth output argument and a sparse input matrix, |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
113 @code{lu} attempts to use a scaling factor @var{R} on the input matrix |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
114 such that |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
115 @code{@var{P} * (@var{R} \ @var{A}) * @var{Q} = @var{L} * @var{U}}. |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
116 This typically leads to a sparser and more stable factorization. |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
117 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
118 An additional input argument @var{thres}, that defines the pivoting |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
119 threshold can be given. @var{thres} can be a scalar, in which case |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
120 it defines the @sc{umfpack} pivoting tolerance for both symmetric and |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
121 unsymmetric cases. If @var{thres} is a 2-element vector, then the first |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
122 element defines the pivoting tolerance for the unsymmetric @sc{umfpack} |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
123 pivoting strategy and the second for the symmetric strategy. By default, |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
124 the values defined by @code{spparms} are used ([0.1, 0.001]). |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
125 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
126 Given the string argument @qcode{"vector"}, @code{lu} returns the values |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
127 of @var{P} and @var{Q} as vector values, such that for full matrix, |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
128 @code{@var{A}(@var{P},:) = @var{L} * @var{U}}, and @code{@var{R}(@var{P},:) |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
129 * @var{A}(:,@var{Q}) = @var{L} * @var{U}}. |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
130 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
131 With two output arguments, returns the permuted forms of the upper and |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
132 lower triangular matrices, such that @code{@var{A} = @var{L} * @var{U}}. |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
133 With one output argument @var{y}, then the matrix returned by the |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
134 @sc{lapack} routines is returned. If the input matrix is sparse then the |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
135 matrix @var{L} is embedded into @var{U} to give a return value similar to |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
136 the full case. For both full and sparse matrices, @code{lu} loses the |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
137 permutation information. |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
138 @seealso{luupdate, ilu, chol, hess, qr, qz, schur, svd} |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
139 @end deftypefn */) |
2928 | 140 { |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
141 int nargin = args.length (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
142 bool issparse = (nargin > 0 && args(0).is_sparse_type ()); |
2928 | 143 |
20909
03e4ddd49396
omit unnecessary nargout checks
John W. Eaton <jwe@octave.org>
parents:
20898
diff
changeset
|
144 if (nargin < 1 || (issparse && nargin > 3) || (! issparse && nargin > 2)) |
20802
8bb38ba1bad6
eliminate return statements after calls to print_usage
John W. Eaton <jwe@octave.org>
parents:
20711
diff
changeset
|
145 print_usage (); |
2928 | 146 |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
147 bool vecout = false; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
148 Matrix thres; |
20892 | 149 int n = 1; |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
150 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
151 while (n < nargin) |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
152 { |
18112
b560bac0fca2
maint: Don't use space between 'args' and '(' when doing indexing.
Rik <rik@octave.org>
parents:
18100
diff
changeset
|
153 if (args(n).is_string ()) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
154 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
155 std::string tmp = args(n++).string_value (); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
156 |
20892 | 157 if (tmp == "vector") |
19743
67f2c76f9f4d
Remove unnecessary checking of error_state after is_string validation.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
158 vecout = true; |
67f2c76f9f4d
Remove unnecessary checking of error_state after is_string validation.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
159 else |
67f2c76f9f4d
Remove unnecessary checking of error_state after is_string validation.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
160 error ("lu: unrecognized string argument"); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
161 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
162 else |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
163 { |
20892 | 164 if (! issparse) |
165 error ("lu: can not define pivoting threshold THRES for full matrices"); | |
21055
5e00ed38a58b
maint: Replace if/error/else paradigm with just if/error.
Rik <rik@octave.org>
parents:
20941
diff
changeset
|
166 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
167 Matrix tmp = args(n++).matrix_value (); |
20892 | 168 if (tmp.numel () == 1) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
169 { |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
170 thres.resize (1,2); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
171 thres(0) = tmp(0); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
172 thres(1) = tmp(0); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
173 } |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
174 else if (tmp.numel () == 2) |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
175 thres = tmp; |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
176 else |
20711
7b608fadc663
Make error messages more specific about the variable and problem encountered.
Rik <rik@octave.org>
parents:
20555
diff
changeset
|
177 error ("lu: THRES must be a 1 or 2-element vector"); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
178 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
179 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
180 |
20892 | 181 octave_value_list retval; |
182 bool scale = (nargout == 5); | |
183 | |
2928 | 184 octave_value arg = args(0); |
185 | |
5275 | 186 octave_idx_type nr = arg.rows (); |
187 octave_idx_type nc = arg.columns (); | |
2928 | 188 |
189 int arg_is_empty = empty_arg ("lu", nr, nc); | |
190 | |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
191 if (issparse) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
192 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
193 if (arg_is_empty < 0) |
20941
a4f5da7c5463
maint: Replace "octave_value_list ()" with "ovl ()".
Rik <rik@octave.org>
parents:
20940
diff
changeset
|
194 return ovl (); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
195 else if (arg_is_empty > 0) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
196 return octave_value_list (5, SparseMatrix ()); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
197 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
198 if (arg.is_real_type ()) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
199 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
200 SparseMatrix m = arg.sparse_matrix_value (); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
201 |
18678
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
18497
diff
changeset
|
202 if (nargout < 4) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
203 { |
20892 | 204 ColumnVector Qinit (nc); |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
205 for (octave_idx_type i = 0; i < nc; i++) |
20892 | 206 Qinit(i) = i; |
21146
ea9c05014809
revamp sparse LU factorization classes
John W. Eaton <jwe@octave.org>
parents:
21139
diff
changeset
|
207 sparse_lu<SparseMatrix> fact (m, Qinit, thres, false, true); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
208 |
18678
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
18497
diff
changeset
|
209 if (nargout < 2) |
19861
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
210 retval(0) = fact.Y (); |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
211 else |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
212 { |
20921
4d3daf7e43f3
eliminate trailing whitespace in source files
John W. Eaton <jwe@octave.org>
parents:
20909
diff
changeset
|
213 retval.resize (nargout == 3 ? 3 : 2); |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
214 retval(1) |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
215 = octave_value ( |
19593
446c46af4b42
strip trailing whitespace from most source files
John W. Eaton <jwe@octave.org>
parents:
18488
diff
changeset
|
216 fact.U () * fact.Pc_mat ().transpose (), |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
217 MatrixType (MatrixType::Permuted_Upper, |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
218 nc, fact.col_perm ())); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
219 |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
220 PermMatrix P = fact.Pr_mat (); |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
221 SparseMatrix L = fact.L (); |
20892 | 222 |
223 if (nargout == 2) | |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
224 retval(0) |
18678
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
18497
diff
changeset
|
225 = octave_value (P.transpose () * L, |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
226 MatrixType (MatrixType::Permuted_Lower, |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
227 nr, fact.row_perm ())); |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
228 else |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
229 { |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
230 retval(0) = L; |
18678
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
18497
diff
changeset
|
231 if (vecout) |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
232 retval(2) = fact.Pr_vec(); |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
233 else |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
234 retval(2) = P; |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
235 } |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
236 } |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
237 } |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
238 else |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
239 { |
20892 | 240 retval.resize (scale ? 5 : 4); |
21146
ea9c05014809
revamp sparse LU factorization classes
John W. Eaton <jwe@octave.org>
parents:
21139
diff
changeset
|
241 sparse_lu<SparseMatrix> fact (m, thres, scale); |
2928 | 242 |
20892 | 243 retval(0) = octave_value (fact.L (), |
244 MatrixType (MatrixType::Lower)); | |
245 retval(1) = octave_value (fact.U (), | |
246 MatrixType (MatrixType::Upper)); | |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
247 |
19861
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
248 if (vecout) |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
249 { |
20892 | 250 retval(2) = fact.Pr_vec (); |
19861
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
251 retval(3) = fact.Pc_vec (); |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
252 } |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
253 else |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
254 { |
20892 | 255 retval(2) = fact.Pr_mat (); |
19861
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
256 retval(3) = fact.Pc_mat (); |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
257 } |
20892 | 258 |
259 if (scale) | |
260 retval(4) = fact.R (); | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
261 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
262 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
263 else if (arg.is_complex_type ()) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
264 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
265 SparseComplexMatrix m = arg.sparse_complex_matrix_value (); |
2928 | 266 |
18678
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
18497
diff
changeset
|
267 if (nargout < 4) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
268 { |
20892 | 269 ColumnVector Qinit (nc); |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
270 for (octave_idx_type i = 0; i < nc; i++) |
20892 | 271 Qinit(i) = i; |
21146
ea9c05014809
revamp sparse LU factorization classes
John W. Eaton <jwe@octave.org>
parents:
21139
diff
changeset
|
272 sparse_lu<SparseComplexMatrix> fact (m, Qinit, thres, false, true); |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
273 |
18678
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
18497
diff
changeset
|
274 if (nargout < 2) |
19861
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
275 retval(0) = fact.Y (); |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
276 else |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
277 { |
20921
4d3daf7e43f3
eliminate trailing whitespace in source files
John W. Eaton <jwe@octave.org>
parents:
20909
diff
changeset
|
278 retval.resize (nargout == 3 ? 3 : 2); |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
279 retval(1) |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
280 = octave_value ( |
19593
446c46af4b42
strip trailing whitespace from most source files
John W. Eaton <jwe@octave.org>
parents:
18488
diff
changeset
|
281 fact.U () * fact.Pc_mat ().transpose (), |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
282 MatrixType (MatrixType::Permuted_Upper, |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
283 nc, fact.col_perm ())); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
284 |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
285 PermMatrix P = fact.Pr_mat (); |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
286 SparseComplexMatrix L = fact.L (); |
20892 | 287 if (nargout == 2) |
19861
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
288 retval(0) |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
289 = octave_value (P.transpose () * L, |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
290 MatrixType (MatrixType::Permuted_Lower, |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
291 nr, fact.row_perm ())); |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
292 else |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
293 { |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
294 retval(0) = L; |
18678
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
18497
diff
changeset
|
295 if (vecout) |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
296 retval(2) = fact.Pr_vec(); |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
297 else |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
298 retval(2) = P; |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
299 } |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
300 } |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
301 } |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
302 else |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
303 { |
20892 | 304 retval.resize (scale ? 5 : 4); |
21146
ea9c05014809
revamp sparse LU factorization classes
John W. Eaton <jwe@octave.org>
parents:
21139
diff
changeset
|
305 sparse_lu<SparseComplexMatrix> fact (m, thres, scale); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
306 |
20892 | 307 retval(0) = octave_value (fact.L (), |
308 MatrixType (MatrixType::Lower)); | |
309 retval(1) = octave_value (fact.U (), | |
310 MatrixType (MatrixType::Upper)); | |
311 | |
312 if (vecout) | |
313 { | |
314 retval(2) = fact.Pr_vec (); | |
315 retval(3) = fact.Pc_vec (); | |
316 } | |
317 else | |
318 { | |
319 retval(2) = fact.Pr_mat (); | |
320 retval(3) = fact.Pc_mat (); | |
321 } | |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
322 |
19861
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
323 if (scale) |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
324 retval(4) = fact.R (); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
325 } |
18427
93c019f00e59
maint: whitespace fixes for libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
18426
diff
changeset
|
326 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
327 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
328 else |
21100
e39e05d90788
Switch gripe_XXX to either err_XXX or warn_XXX naming scheme.
Rik <rik@octave.org>
parents:
21055
diff
changeset
|
329 err_wrong_type_arg ("lu", arg); |
2928 | 330 } |
331 else | |
332 { | |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
333 if (arg_is_empty < 0) |
20941
a4f5da7c5463
maint: Replace "octave_value_list ()" with "ovl ()".
Rik <rik@octave.org>
parents:
20940
diff
changeset
|
334 return ovl (); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
335 else if (arg_is_empty > 0) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
336 return octave_value_list (3, Matrix ()); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
337 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
338 if (arg.is_real_type ()) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
339 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
340 if (arg.is_single_type ()) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
341 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
342 FloatMatrix m = arg.float_matrix_value (); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
343 |
21271
7e67c7f82fc1
better use of templates for lu factorization classes
John W. Eaton <jwe@octave.org>
parents:
21200
diff
changeset
|
344 lu<FloatMatrix> fact (m); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
345 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
346 switch (nargout) |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
347 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
348 case 0: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
349 case 1: |
20892 | 350 retval = ovl (fact.Y ()); |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
351 break; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
352 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
353 case 2: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
354 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
355 PermMatrix P = fact.P (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
356 FloatMatrix L = P.transpose () * fact.L (); |
20892 | 357 retval = ovl (L, get_lu_u (fact)); |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
358 } |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
359 break; |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
360 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
361 case 3: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
362 default: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
363 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
364 if (vecout) |
20892 | 365 retval = ovl (get_lu_l (fact), get_lu_u (fact), |
366 fact.P_vec ()); | |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
367 else |
20892 | 368 retval = ovl (get_lu_l (fact), get_lu_u (fact), |
369 fact.P ()); | |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
370 } |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
371 break; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
372 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
373 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
374 else |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
375 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
376 Matrix m = arg.matrix_value (); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
377 |
21271
7e67c7f82fc1
better use of templates for lu factorization classes
John W. Eaton <jwe@octave.org>
parents:
21200
diff
changeset
|
378 lu<Matrix> fact (m); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
379 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
380 switch (nargout) |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
381 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
382 case 0: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
383 case 1: |
20892 | 384 retval = ovl (fact.Y ()); |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
385 break; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
386 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
387 case 2: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
388 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
389 PermMatrix P = fact.P (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
390 Matrix L = P.transpose () * fact.L (); |
20892 | 391 retval = ovl (L, get_lu_u (fact)); |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
392 } |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
393 break; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
394 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
395 case 3: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
396 default: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
397 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
398 if (vecout) |
20892 | 399 retval = ovl (get_lu_l (fact), get_lu_u (fact), |
400 fact.P_vec ()); | |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
401 else |
20892 | 402 retval = ovl (get_lu_l (fact), get_lu_u (fact), |
403 fact.P ()); | |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
404 } |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
405 break; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
406 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
407 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
408 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
409 else if (arg.is_complex_type ()) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
410 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
411 if (arg.is_single_type ()) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
412 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
413 FloatComplexMatrix m = arg.float_complex_matrix_value (); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
414 |
21271
7e67c7f82fc1
better use of templates for lu factorization classes
John W. Eaton <jwe@octave.org>
parents:
21200
diff
changeset
|
415 lu<FloatComplexMatrix> fact (m); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
416 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
417 switch (nargout) |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
418 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
419 case 0: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
420 case 1: |
20892 | 421 retval = ovl (fact.Y ()); |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
422 break; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
423 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
424 case 2: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
425 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
426 PermMatrix P = fact.P (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
427 FloatComplexMatrix L = P.transpose () * fact.L (); |
20892 | 428 retval = ovl (L, get_lu_u (fact)); |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
429 } |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
430 break; |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
431 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
432 case 3: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
433 default: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
434 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
435 if (vecout) |
20892 | 436 retval = ovl (get_lu_l (fact), get_lu_u (fact), |
437 fact.P_vec ()); | |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
438 else |
20892 | 439 retval = ovl (get_lu_l (fact), get_lu_u (fact), |
440 fact.P ()); | |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
441 } |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
442 break; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
443 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
444 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
445 else |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
446 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
447 ComplexMatrix m = arg.complex_matrix_value (); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
448 |
21271
7e67c7f82fc1
better use of templates for lu factorization classes
John W. Eaton <jwe@octave.org>
parents:
21200
diff
changeset
|
449 lu<ComplexMatrix> fact (m); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
450 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
451 switch (nargout) |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
452 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
453 case 0: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
454 case 1: |
20892 | 455 retval = ovl (fact.Y ()); |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
456 break; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
457 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
458 case 2: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
459 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
460 PermMatrix P = fact.P (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
461 ComplexMatrix L = P.transpose () * fact.L (); |
20892 | 462 retval = ovl (L, get_lu_u (fact)); |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
463 } |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
464 break; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
465 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
466 case 3: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
467 default: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
468 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
469 if (vecout) |
20892 | 470 retval = ovl (get_lu_l (fact), get_lu_u (fact), |
471 fact.P_vec ()); | |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
472 else |
20892 | 473 retval = ovl (get_lu_l (fact), get_lu_u (fact), |
474 fact.P ()); | |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
475 } |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
476 break; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
477 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
478 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
479 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
480 else |
21100
e39e05d90788
Switch gripe_XXX to either err_XXX or warn_XXX naming scheme.
Rik <rik@octave.org>
parents:
21055
diff
changeset
|
481 err_wrong_type_arg ("lu", arg); |
2928 | 482 } |
483 | |
484 return retval; | |
485 } | |
486 | |
487 /* | |
21317
a4faec57f4c8
maint: remove semicolon after %!assert tests to follow Octave conventions.
Rik <rik@octave.org>
parents:
21301
diff
changeset
|
488 %!assert(lu ([1, 2; 3, 4]), [3, 4; 1/3, 2/3], eps) |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
489 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
490 %!test |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
491 %! [l, u] = lu ([1, 2; 3, 4]); |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
492 %! assert (l, [1/3, 1; 1, 0], sqrt (eps)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
493 %! assert (u, [3, 4; 0, 2/3], sqrt (eps)); |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
494 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
495 %!test |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
496 %! [l, u, p] = lu ([1, 2; 3, 4]); |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
497 %! assert (l, [1, 0; 1/3, 1], sqrt (eps)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
498 %! assert (u, [3, 4; 0, 2/3], sqrt (eps)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
499 %! assert (p(:,:), [0, 1; 1, 0], sqrt (eps)); |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
500 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
501 %!test |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
502 %! [l, u, p] = lu ([1, 2; 3, 4], "vector"); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
503 %! assert (l, [1, 0; 1/3, 1], sqrt (eps)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
504 %! assert (u, [3, 4; 0, 2/3], sqrt (eps)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
505 %! assert (p, [2;1], sqrt (eps)); |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
506 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
507 %!test |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
508 %! [l, u, p] = lu ([1, 2; 3, 4; 5, 6]); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
509 %! assert (l, [1, 0; 1/5, 1; 3/5, 1/2], sqrt (eps)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
510 %! assert (u, [5, 6; 0, 4/5], sqrt (eps)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
511 %! assert (p(:,:), [0, 0, 1; 1, 0, 0; 0 1 0], sqrt (eps)); |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
512 |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
513 %!assert (lu (single ([1, 2; 3, 4])), single ([3, 4; 1/3, 2/3]), eps ("single")) |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
514 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
515 %!test |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
516 %! [l, u] = lu (single ([1, 2; 3, 4])); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
517 %! assert (l, single ([1/3, 1; 1, 0]), sqrt (eps ("single"))); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
518 %! assert (u, single ([3, 4; 0, 2/3]), sqrt (eps ("single"))); |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
519 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
520 %!test |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
521 %! [l, u, p] = lu (single ([1, 2; 3, 4])); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
522 %! assert (l, single ([1, 0; 1/3, 1]), sqrt (eps ("single"))); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
523 %! assert (u, single ([3, 4; 0, 2/3]), sqrt (eps ("single"))); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
524 %! assert (p(:,:), single ([0, 1; 1, 0]), sqrt (eps ("single"))); |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
525 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
526 %!test |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
527 %! [l, u, p] = lu (single ([1, 2; 3, 4]), "vector"); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
528 %! assert (l, single ([1, 0; 1/3, 1]), sqrt (eps ("single"))); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
529 %! assert (u, single ([3, 4; 0, 2/3]), sqrt (eps ("single"))); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
530 %! assert (p, single ([2;1]), sqrt (eps ("single"))); |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
531 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
532 %!test |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
533 %! [l u p] = lu (single ([1, 2; 3, 4; 5, 6])); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
534 %! assert (l, single ([1, 0; 1/5, 1; 3/5, 1/2]), sqrt (eps ("single"))); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
535 %! assert (u, single ([5, 6; 0, 4/5]), sqrt (eps ("single"))); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
536 %! assert (p(:,:), single ([0, 0, 1; 1, 0, 0; 0 1 0]), sqrt (eps ("single"))); |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
537 |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
538 %!error lu () |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
539 %!error <can not define pivoting threshold> lu ([1, 2; 3, 4], 2) |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
540 |
18488
4c2465444a96
Use %!testif HAVE_UMFPACK for sparse lu test added in cset 2a45b6b87bee
Mike Miller <mtmiller@ieee.org>
parents:
18427
diff
changeset
|
541 %!testif HAVE_UMFPACK |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
542 %! Bi = [1 2 3 4 5 2 3 6 7 8 4 5 7 8 9]; |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
543 %! Bj = [1 3 4 5 6 7 8 9 11 12 13 14 15 16 17]; |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
544 %! Bv = [1 1 1 1 1 1 -1 1 1 1 1 -1 1 -1 1]; |
18427
93c019f00e59
maint: whitespace fixes for libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
18426
diff
changeset
|
545 %! B = sparse (Bi, Bj, Bv); |
93c019f00e59
maint: whitespace fixes for libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
18426
diff
changeset
|
546 %! [L, U] = lu (B); |
18426
807e21ebddd5
maint: Style fixes for tests in libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
18425
diff
changeset
|
547 %! assert (L*U, B); |
18427
93c019f00e59
maint: whitespace fixes for libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
18426
diff
changeset
|
548 %! [L, U, P] = lu(B); |
18426
807e21ebddd5
maint: Style fixes for tests in libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
18425
diff
changeset
|
549 %! assert (P'*L*U, B); |
18427
93c019f00e59
maint: whitespace fixes for libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
18426
diff
changeset
|
550 %! [L, U, P, Q] = lu (B); |
18426
807e21ebddd5
maint: Style fixes for tests in libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
18425
diff
changeset
|
551 %! assert (P'*L*U*Q', B); |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
552 |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
553 */ |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
554 |
9708 | 555 static |
556 bool check_lu_dims (const octave_value& l, const octave_value& u, | |
557 const octave_value& p) | |
558 { | |
18100
6a71e5030df5
Follow coding convention of defining and initializing only 1 variable per line in liboctinterp.
Rik <rik@octave.org>
parents:
17787
diff
changeset
|
559 octave_idx_type m = l.rows (); |
6a71e5030df5
Follow coding convention of defining and initializing only 1 variable per line in liboctinterp.
Rik <rik@octave.org>
parents:
17787
diff
changeset
|
560 octave_idx_type k = u.rows (); |
6a71e5030df5
Follow coding convention of defining and initializing only 1 variable per line in liboctinterp.
Rik <rik@octave.org>
parents:
17787
diff
changeset
|
561 octave_idx_type n = u.columns (); |
20892 | 562 |
9708 | 563 return ((l.ndims () == 2 && u.ndims () == 2 && k == l.columns ()) |
19864
17d647821d61
maint: More cleanup of C++ code to follow Octave coding conventions.
John W. Eaton <jwe@octave.org>
parents:
19861
diff
changeset
|
564 && k == std::min (m, n) |
17d647821d61
maint: More cleanup of C++ code to follow Octave coding conventions.
John W. Eaton <jwe@octave.org>
parents:
19861
diff
changeset
|
565 && (p.is_undefined () || p.rows () == m)); |
9708 | 566 } |
567 | |
15039
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14854
diff
changeset
|
568 DEFUN (luupdate, args, , |
21966
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
569 doc: /* -*- texinfo -*- |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
570 @deftypefn {} {[@var{L}, @var{U}] =} luupdate (@var{L}, @var{U}, @var{x}, @var{y}) |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
571 @deftypefnx {} {[@var{L}, @var{U}, @var{P}] =} luupdate (@var{L}, @var{U}, @var{P}, @var{x}, @var{y}) |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
572 Given an LU@tie{}factorization of a real or complex matrix |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
573 @w{@var{A} = @var{L}*@var{U}}, @var{L}@tie{}lower unit trapezoidal and |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
574 @var{U}@tie{}upper trapezoidal, return the LU@tie{}factorization |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
575 of @w{@var{A} + @var{x}*@var{y}.'}, where @var{x} and @var{y} are |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
576 column vectors (rank-1 update) or matrices with equal number of columns |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
577 (rank-k update). |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
578 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
579 Optionally, row-pivoted updating can be used by supplying a row permutation |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
580 (pivoting) matrix @var{P}; in that case, an updated permutation matrix is |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
581 returned. Note that if @var{L}, @var{U}, @var{P} is a pivoted |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
582 LU@tie{}factorization as obtained by @code{lu}: |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
583 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
584 @example |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
585 [@var{L}, @var{U}, @var{P}] = lu (@var{A}); |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
586 @end example |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
587 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
588 @noindent |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
589 then a factorization of @tcode{@var{A}+@var{x}*@var{y}.'} can be obtained |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
590 either as |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
591 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
592 @example |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
593 [@var{L1}, @var{U1}] = lu (@var{L}, @var{U}, @var{P}*@var{x}, @var{y}) |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
594 @end example |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
595 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
596 @noindent |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
597 or |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
598 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
599 @example |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
600 [@var{L1}, @var{U1}, @var{P1}] = lu (@var{L}, @var{U}, @var{P}, @var{x}, @var{y}) |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
601 @end example |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
602 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
603 The first form uses the unpivoted algorithm, which is faster, but less |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
604 stable. The second form uses a slower pivoted algorithm, which is more |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
605 stable. |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
606 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
607 The matrix case is done as a sequence of rank-1 updates; thus, for large |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
608 enough k, it will be both faster and more accurate to recompute the |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
609 factorization from scratch. |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
610 @seealso{lu, cholupdate, qrupdate} |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
611 @end deftypefn */) |
9708 | 612 { |
20816
b16bcd7a2a33
Use int rather than octave_idx_type for nargin data type.
Rik <rik@octave.org>
parents:
20812
diff
changeset
|
613 int nargin = args.length (); |
9708 | 614 |
615 if (nargin != 4 && nargin != 5) | |
20802
8bb38ba1bad6
eliminate return statements after calls to print_usage
John W. Eaton <jwe@octave.org>
parents:
20711
diff
changeset
|
616 print_usage (); |
9708 | 617 |
20812
d9ca869ca124
maint: Clean-up more instances of print_usage().
Rik <rik@octave.org>
parents:
20802
diff
changeset
|
618 bool pivoted = (nargin == 5); |
d9ca869ca124
maint: Clean-up more instances of print_usage().
Rik <rik@octave.org>
parents:
20802
diff
changeset
|
619 |
9708 | 620 octave_value argl = args(0); |
621 octave_value argu = args(1); | |
622 octave_value argp = pivoted ? args(2) : octave_value (); | |
623 octave_value argx = args(2 + pivoted); | |
624 octave_value argy = args(3 + pivoted); | |
625 | |
20892 | 626 if (! (argl.is_numeric_type () && argu.is_numeric_type () |
627 && argx.is_numeric_type () && argy.is_numeric_type () | |
628 && (! pivoted || argp.is_perm_matrix ()))) | |
629 error ("luupdate: L, U, X, and Y must be numeric"); | |
9708 | 630 |
20892 | 631 if (! check_lu_dims (argl, argu, argp)) |
632 error ("luupdate: dimension mismatch"); | |
9708 | 633 |
20892 | 634 PermMatrix P = (pivoted |
635 ? argp.perm_matrix_value () | |
636 : PermMatrix::eye (argl.rows ())); | |
9708 | 637 |
20892 | 638 if (argl.is_real_type () && argu.is_real_type () |
639 && argx.is_real_type () && argy.is_real_type ()) | |
640 { | |
641 // all real case | |
642 if (argl.is_single_type () || argu.is_single_type () | |
643 || argx.is_single_type () || argy.is_single_type ()) | |
644 { | |
645 FloatMatrix L = argl.float_matrix_value (); | |
646 FloatMatrix U = argu.float_matrix_value (); | |
647 FloatMatrix x = argx.float_matrix_value (); | |
648 FloatMatrix y = argy.float_matrix_value (); | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
649 |
21271
7e67c7f82fc1
better use of templates for lu factorization classes
John W. Eaton <jwe@octave.org>
parents:
21200
diff
changeset
|
650 lu<FloatMatrix> fact (L, U, P); |
20892 | 651 if (pivoted) |
652 fact.update_piv (x, y); | |
653 else | |
654 fact.update (x, y); | |
9708 | 655 |
20892 | 656 if (pivoted) |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20892
diff
changeset
|
657 return ovl (get_lu_l (fact), get_lu_u (fact), fact.P ()); |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20892
diff
changeset
|
658 else |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20892
diff
changeset
|
659 return ovl (get_lu_l (fact), get_lu_u (fact)); |
9708 | 660 } |
661 else | |
20892 | 662 { |
663 Matrix L = argl.matrix_value (); | |
664 Matrix U = argu.matrix_value (); | |
665 Matrix x = argx.matrix_value (); | |
666 Matrix y = argy.matrix_value (); | |
667 | |
21271
7e67c7f82fc1
better use of templates for lu factorization classes
John W. Eaton <jwe@octave.org>
parents:
21200
diff
changeset
|
668 lu<Matrix> fact (L, U, P); |
20892 | 669 if (pivoted) |
670 fact.update_piv (x, y); | |
671 else | |
672 fact.update (x, y); | |
673 | |
674 if (pivoted) | |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20892
diff
changeset
|
675 return ovl (get_lu_l (fact), get_lu_u (fact), fact.P ()); |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20892
diff
changeset
|
676 else |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20892
diff
changeset
|
677 return ovl (get_lu_l (fact), get_lu_u (fact)); |
20892 | 678 } |
9708 | 679 } |
680 else | |
20892 | 681 { |
682 // complex case | |
683 if (argl.is_single_type () || argu.is_single_type () | |
684 || argx.is_single_type () || argy.is_single_type ()) | |
685 { | |
686 FloatComplexMatrix L = argl.float_complex_matrix_value (); | |
687 FloatComplexMatrix U = argu.float_complex_matrix_value (); | |
688 FloatComplexMatrix x = argx.float_complex_matrix_value (); | |
689 FloatComplexMatrix y = argy.float_complex_matrix_value (); | |
690 | |
21271
7e67c7f82fc1
better use of templates for lu factorization classes
John W. Eaton <jwe@octave.org>
parents:
21200
diff
changeset
|
691 lu<FloatComplexMatrix> fact (L, U, P); |
20892 | 692 if (pivoted) |
693 fact.update_piv (x, y); | |
694 else | |
695 fact.update (x, y); | |
696 | |
697 if (pivoted) | |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20892
diff
changeset
|
698 return ovl (get_lu_l (fact), get_lu_u (fact), fact.P ()); |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20892
diff
changeset
|
699 else |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20892
diff
changeset
|
700 return ovl (get_lu_l (fact), get_lu_u (fact)); |
20892 | 701 } |
702 else | |
703 { | |
704 ComplexMatrix L = argl.complex_matrix_value (); | |
705 ComplexMatrix U = argu.complex_matrix_value (); | |
706 ComplexMatrix x = argx.complex_matrix_value (); | |
707 ComplexMatrix y = argy.complex_matrix_value (); | |
708 | |
21271
7e67c7f82fc1
better use of templates for lu factorization classes
John W. Eaton <jwe@octave.org>
parents:
21200
diff
changeset
|
709 lu<ComplexMatrix> fact (L, U, P); |
20892 | 710 if (pivoted) |
711 fact.update_piv (x, y); | |
712 else | |
713 fact.update (x, y); | |
714 | |
715 if (pivoted) | |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20892
diff
changeset
|
716 return ovl (get_lu_l (fact), get_lu_u (fact), fact.P ()); |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20892
diff
changeset
|
717 else |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20892
diff
changeset
|
718 return ovl (get_lu_l (fact), get_lu_u (fact)); |
20892 | 719 } |
720 } | |
9708 | 721 } |
722 | |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
723 /* |
10080
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
724 %!shared A, u, v, Ac, uc, vc |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
725 %! A = [0.091364 0.613038 0.999083; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
726 %! 0.594638 0.425302 0.603537; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
727 %! 0.383594 0.291238 0.085574; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
728 %! 0.265712 0.268003 0.238409; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
729 %! 0.669966 0.743851 0.445057 ]; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
730 %! |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
731 %! u = [0.85082; |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
732 %! 0.76426; |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
733 %! 0.42883; |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
734 %! 0.53010; |
10080
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
735 %! 0.80683 ]; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
736 %! |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
737 %! v = [0.98810; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
738 %! 0.24295; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
739 %! 0.43167 ]; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
740 %! |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
741 %! Ac = [0.620405 + 0.956953i 0.480013 + 0.048806i 0.402627 + 0.338171i; |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
742 %! 0.589077 + 0.658457i 0.013205 + 0.279323i 0.229284 + 0.721929i; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
743 %! 0.092758 + 0.345687i 0.928679 + 0.241052i 0.764536 + 0.832406i; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
744 %! 0.912098 + 0.721024i 0.049018 + 0.269452i 0.730029 + 0.796517i; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
745 %! 0.112849 + 0.603871i 0.486352 + 0.142337i 0.355646 + 0.151496i ]; |
10080
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
746 %! |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
747 %! uc = [0.20351 + 0.05401i; |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
748 %! 0.13141 + 0.43708i; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
749 %! 0.29808 + 0.08789i; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
750 %! 0.69821 + 0.38844i; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
751 %! 0.74871 + 0.25821i ]; |
10080
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
752 %! |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
753 %! vc = [0.85839 + 0.29468i; |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
754 %! 0.20820 + 0.93090i; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
755 %! 0.86184 + 0.34689i ]; |
10080
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
756 %! |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
757 |
10094
ab1101011a6d
lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents:
10080
diff
changeset
|
758 %!testif HAVE_QRUPDATE_LUU |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
759 %! [L,U,P] = lu (A); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
760 %! [L,U] = luupdate (L,U,P*u,v); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
761 %! assert (norm (vec (tril (L)-L), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
762 %! assert (norm (vec (triu (U)-U), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
763 %! assert (norm (vec (P'*L*U - A - u*v.'), Inf) < norm (A)*1e1*eps); |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
764 %! |
10094
ab1101011a6d
lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents:
10080
diff
changeset
|
765 %!testif HAVE_QRUPDATE_LUU |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
766 %! [L,U,P] = lu (Ac); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
767 %! [L,U] = luupdate (L,U,P*uc,vc); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
768 %! assert (norm (vec (tril (L)-L), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
769 %! assert (norm (vec (triu (U)-U), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
770 %! assert (norm (vec (P'*L*U - Ac - uc*vc.'), Inf) < norm (Ac)*1e1*eps); |
10080
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
771 |
10094
ab1101011a6d
lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents:
10080
diff
changeset
|
772 %!testif HAVE_QRUPDATE_LUU |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
773 %! [L,U,P] = lu (single (A)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
774 %! [L,U] = luupdate (L,U,P*single (u), single (v)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
775 %! assert (norm (vec (tril (L)-L), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
776 %! assert (norm (vec (triu (U)-U), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
777 %! assert (norm (vec (P'*L*U - single (A) - single (u)*single (v).'), Inf) < norm (single (A))*1e1*eps ("single")); |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
778 %! |
10094
ab1101011a6d
lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents:
10080
diff
changeset
|
779 %!testif HAVE_QRUPDATE_LUU |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
780 %! [L,U,P] = lu (single (Ac)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
781 %! [L,U] = luupdate (L,U,P*single (uc),single (vc)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
782 %! assert (norm (vec (tril (L)-L), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
783 %! assert (norm (vec (triu (U)-U), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
784 %! assert (norm (vec (P'*L*U - single (Ac) - single (uc)*single (vc).'), Inf) < norm (single (Ac))*1e1*eps ("single")); |
10080
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
785 |
10094
ab1101011a6d
lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents:
10080
diff
changeset
|
786 %!testif HAVE_QRUPDATE_LUU |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
787 %! [L,U,P] = lu (A); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
788 %! [L,U,P] = luupdate (L,U,P,u,v); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
789 %! assert (norm (vec (tril (L)-L), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
790 %! assert (norm (vec (triu (U)-U), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
791 %! assert (norm (vec (P'*L*U - A - u*v.'), Inf) < norm (A)*1e1*eps); |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
792 %! |
10094
ab1101011a6d
lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents:
10080
diff
changeset
|
793 %!testif HAVE_QRUPDATE_LUU |
18849
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18678
diff
changeset
|
794 %! [L,U,P] = lu (A); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18678
diff
changeset
|
795 %! [~,ordcols] = max (P,[],1); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18678
diff
changeset
|
796 %! [~,ordrows] = max (P,[],2); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18678
diff
changeset
|
797 %! P1 = eye (size(P))(:,ordcols); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18678
diff
changeset
|
798 %! P2 = eye (size(P))(ordrows,:); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18678
diff
changeset
|
799 %! assert(P1 == P); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18678
diff
changeset
|
800 %! assert(P2 == P); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18678
diff
changeset
|
801 %! [L,U,P] = luupdate (L,U,P,u,v); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18678
diff
changeset
|
802 %! [L,U,P1] = luupdate (L,U,P1,u,v); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18678
diff
changeset
|
803 %! [L,U,P2] = luupdate (L,U,P2,u,v); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18678
diff
changeset
|
804 %! assert(P1 == P); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18678
diff
changeset
|
805 %! assert(P2 == P); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18678
diff
changeset
|
806 %! |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18678
diff
changeset
|
807 %!testif HAVE_QRUPDATE_LUU |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
808 %! [L,U,P] = lu (Ac); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
809 %! [L,U,P] = luupdate (L,U,P,uc,vc); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
810 %! assert (norm (vec (tril (L)-L), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
811 %! assert (norm (vec (triu (U)-U), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
812 %! assert (norm (vec (P'*L*U - Ac - uc*vc.'), Inf) < norm (Ac)*1e1*eps); |
10080
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
813 |
10094
ab1101011a6d
lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents:
10080
diff
changeset
|
814 %!testif HAVE_QRUPDATE_LUU |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
815 %! [L,U,P] = lu (single (A)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
816 %! [L,U,P] = luupdate (L,U,P,single (u),single (v)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
817 %! assert (norm (vec (tril (L)-L), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
818 %! assert (norm (vec (triu (U)-U), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
819 %! assert (norm (vec (P'*L*U - single (A) - single (u)*single (v).'), Inf) < norm (single (A))*1e1*eps ("single")); |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
820 %! |
10094
ab1101011a6d
lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents:
10080
diff
changeset
|
821 %!testif HAVE_QRUPDATE_LUU |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
822 %! [L,U,P] = lu (single (Ac)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
823 %! [L,U,P] = luupdate (L,U,P,single (uc),single (vc)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
824 %! assert (norm (vec (tril (L)-L), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
825 %! assert (norm (vec (triu (U)-U), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
826 %! assert (norm (vec (P'*L*U - single (Ac) - single (uc)*single (vc).'), Inf) < norm (single (Ac))*1e1*eps ("single")); |
10080
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
827 */ |