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