Mercurial > hg > octave-nkf
annotate libinterp/corefcn/lu.cc @ 20830:b65888ec820e draft default tip gccjit
dmalcom gcc jit import
author | Stefan Mahr <dac922@gmx.de> |
---|---|
date | Fri, 27 Feb 2015 16:59:36 +0100 |
parents | f90c8372b7ba |
children |
rev | line source |
---|---|
2928 | 1 /* |
2 | |
19898
4197fc428c7d
maint: Update copyright notices for 2015.
John W. Eaton <jwe@octave.org>
parents:
19793
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\ |
15039
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14854
diff
changeset
|
66 @deftypefn {Built-in Function} {[@var{L}, @var{U}] =} lu (@var{A})\n\ |
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14854
diff
changeset
|
67 @deftypefnx {Built-in Function} {[@var{L}, @var{U}, @var{P}] =} lu (@var{A})\n\ |
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14854
diff
changeset
|
68 @deftypefnx {Built-in Function} {[@var{L}, @var{U}, @var{P}, @var{Q}] =} lu (@var{S})\n\ |
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14854
diff
changeset
|
69 @deftypefnx {Built-in Function} {[@var{L}, @var{U}, @var{P}, @var{Q}, @var{R}] =} lu (@var{S})\n\ |
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14854
diff
changeset
|
70 @deftypefnx {Built-in Function} {[@dots{}] =} lu (@var{S}, @var{thres})\n\ |
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14854
diff
changeset
|
71 @deftypefnx {Built-in Function} {@var{y} =} lu (@dots{})\n\ |
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14854
diff
changeset
|
72 @deftypefnx {Built-in Function} {[@dots{}] =} lu (@dots{}, \"vector\")\n\ |
3372 | 73 @cindex LU decomposition\n\ |
20382
4f45eaf83908
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
20071
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:
20071
diff
changeset
|
75 \n\ |
4f45eaf83908
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
20071
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:
20071
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:
20071
diff
changeset
|
78 \n\ |
4f45eaf83908
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
20071
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\ |
19247
38937efbee21
Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents:
19039
diff
changeset
|
142 @seealso{luupdate, ilu, chol, hess, qr, qz, schur, svd}\n\ |
3372 | 143 @end deftypefn") |
2928 | 144 { |
145 octave_value_list retval; | |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
146 int nargin = args.length (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
147 bool issparse = (nargin > 0 && args(0).is_sparse_type ()); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
148 bool scale = (nargout == 5); |
2928 | 149 |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
150 if (nargin < 1 || (issparse && (nargin > 3 || nargout > 5)) |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
151 || (!issparse && (nargin > 2 || nargout > 3))) |
2928 | 152 { |
5823 | 153 print_usage (); |
2928 | 154 return retval; |
155 } | |
156 | |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
157 bool vecout = false; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
158 Matrix thres; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
159 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
160 int n = 1; |
20762
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
161 while (n < nargin) |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
162 { |
18111
b560bac0fca2
maint: Don't use space between 'args' and '(' when doing indexing.
Rik <rik@octave.org>
parents:
18099
diff
changeset
|
163 if (args(n).is_string ()) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
164 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
165 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
|
166 |
19948
67f2c76f9f4d
Remove unnecessary checking of error_state after is_string validation.
Rik <rik@octave.org>
parents:
19898
diff
changeset
|
167 if (tmp.compare ("vector") == 0) |
67f2c76f9f4d
Remove unnecessary checking of error_state after is_string validation.
Rik <rik@octave.org>
parents:
19898
diff
changeset
|
168 vecout = true; |
67f2c76f9f4d
Remove unnecessary checking of error_state after is_string validation.
Rik <rik@octave.org>
parents:
19898
diff
changeset
|
169 else |
67f2c76f9f4d
Remove unnecessary checking of error_state after is_string validation.
Rik <rik@octave.org>
parents:
19898
diff
changeset
|
170 error ("lu: unrecognized string argument"); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
171 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
172 else |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
173 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
174 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
|
175 |
20762
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
176 if (!issparse) |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
177 error ("lu: can not define pivoting threshold THRES for full matrices"); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
178 else if (tmp.numel () == 1) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
179 { |
20762
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
180 thres.resize (1,2); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
181 thres(0) = tmp(0); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
182 thres(1) = tmp(0); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
183 } |
20762
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
184 else if (tmp.numel () == 2) |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
185 thres = tmp; |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
186 else |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
187 error ("lu: expecting 2-element vector for THRES"); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
188 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
189 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
190 |
2928 | 191 octave_value arg = args(0); |
192 | |
5275 | 193 octave_idx_type nr = arg.rows (); |
194 octave_idx_type nc = arg.columns (); | |
2928 | 195 |
196 int arg_is_empty = empty_arg ("lu", nr, nc); | |
197 | |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
198 if (issparse) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
199 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
200 if (arg_is_empty < 0) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
201 return retval; |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
202 else if (arg_is_empty > 0) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
203 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
|
204 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
205 if (arg.is_real_type ()) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
206 { |
18493
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
207 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
208 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
|
209 |
18833
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
18565
diff
changeset
|
210 if (nargout < 4) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
211 { |
18493
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
212 |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
213 ColumnVector Qinit; |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
214 Qinit.resize (nc); |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
215 for (octave_idx_type i = 0; i < nc; i++) |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
216 Qinit (i) = i; |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
217 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
|
218 |
18833
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
18565
diff
changeset
|
219 if (nargout < 2) |
20068
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
220 retval(0) = fact.Y (); |
18493
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
221 else |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
222 { |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
223 |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
224 retval(1) |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
225 = octave_value ( |
19790
446c46af4b42
strip trailing whitespace from most source files
John W. Eaton <jwe@octave.org>
parents:
18556
diff
changeset
|
226 fact.U () * fact.Pc_mat ().transpose (), |
18493
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
227 MatrixType (MatrixType::Permuted_Upper, |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
228 nc, fact.col_perm ())); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
229 |
18493
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
230 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
|
231 SparseMatrix L = fact.L (); |
18833
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
18565
diff
changeset
|
232 if (nargout < 3) |
18493
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
233 retval(0) |
18833
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
18565
diff
changeset
|
234 = octave_value (P.transpose () * L, |
18493
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
235 MatrixType (MatrixType::Permuted_Lower, |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
236 nr, fact.row_perm ())); |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
237 else |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
238 { |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
239 retval(0) = L; |
18833
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
18565
diff
changeset
|
240 if (vecout) |
18493
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
241 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
|
242 else |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
243 retval(2) = P; |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
244 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
245 |
18493
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
246 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
247 |
18493
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
248 } |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
249 else |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
250 { |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
251 |
20068
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
252 SparseLU fact (m, thres, scale); |
2928 | 253 |
20068
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
254 if (scale) |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
255 retval(4) = fact.R (); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
256 |
20068
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
257 if (vecout) |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
258 { |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
259 retval(3) = fact.Pc_vec (); |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
260 retval(2) = fact.Pr_vec (); |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
261 } |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
262 else |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
263 { |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
264 retval(3) = fact.Pc_mat (); |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
265 retval(2) = fact.Pr_mat (); |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
266 } |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
267 retval(1) = octave_value (fact.U (), |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
268 MatrixType (MatrixType::Upper)); |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
269 retval(0) = octave_value (fact.L (), |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
270 MatrixType (MatrixType::Lower)); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
271 } |
18493
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
272 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
273 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
274 else if (arg.is_complex_type ()) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
275 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
276 SparseComplexMatrix m = arg.sparse_complex_matrix_value (); |
2928 | 277 |
18833
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
18565
diff
changeset
|
278 if (nargout < 4) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
279 { |
18493
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
280 |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
281 ColumnVector Qinit; |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
282 Qinit.resize (nc); |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
283 for (octave_idx_type i = 0; i < nc; i++) |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
284 Qinit (i) = i; |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
285 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
|
286 |
18833
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
18565
diff
changeset
|
287 if (nargout < 2) |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
288 |
20068
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
289 retval(0) = fact.Y (); |
18493
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
290 |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
291 else |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
292 { |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
293 |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
294 retval(1) |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
295 = octave_value ( |
19790
446c46af4b42
strip trailing whitespace from most source files
John W. Eaton <jwe@octave.org>
parents:
18556
diff
changeset
|
296 fact.U () * fact.Pc_mat ().transpose (), |
18493
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
297 MatrixType (MatrixType::Permuted_Upper, |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
298 nc, fact.col_perm ())); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
299 |
18493
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
300 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
|
301 SparseComplexMatrix L = fact.L (); |
18833
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
18565
diff
changeset
|
302 if (nargout < 3) |
20068
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
303 retval(0) |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
304 = octave_value (P.transpose () * L, |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
305 MatrixType (MatrixType::Permuted_Lower, |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
306 nr, fact.row_perm ())); |
18493
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 { |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
309 retval(0) = L; |
18833
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
18565
diff
changeset
|
310 if (vecout) |
18493
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
311 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
|
312 else |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
313 retval(2) = P; |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
314 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
315 |
18493
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
316 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
317 |
18493
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
318 } |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
319 else |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
320 { |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
321 |
20068
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
322 SparseComplexLU fact (m, thres, scale); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
323 |
20068
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
324 if (scale) |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
325 retval(4) = fact.R (); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
326 |
20068
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
327 if (vecout) |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
328 { |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
329 retval(3) = fact.Pc_vec (); |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
330 retval(2) = fact.Pr_vec (); |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
331 } |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
332 else |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
333 { |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
334 retval(3) = fact.Pc_mat (); |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
335 retval(2) = fact.Pr_mat (); |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
336 } |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
337 retval(1) = octave_value (fact.U (), |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
338 MatrixType (MatrixType::Upper)); |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
339 retval(0) = octave_value (fact.L (), |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19948
diff
changeset
|
340 MatrixType (MatrixType::Lower)); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
341 } |
18495
93c019f00e59
maint: whitespace fixes for libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
18494
diff
changeset
|
342 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
343 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
344 else |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
345 gripe_wrong_type_arg ("lu", arg); |
2928 | 346 } |
347 else | |
348 { | |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
349 if (arg_is_empty < 0) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
350 return retval; |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
351 else if (arg_is_empty > 0) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
352 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
|
353 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
354 if (arg.is_real_type ()) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
355 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
356 if (arg.is_single_type ()) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
357 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
358 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
|
359 |
20762
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
360 FloatLU fact (m); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
361 |
20762
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
362 switch (nargout) |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
363 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
364 case 0: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
365 case 1: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
366 retval(0) = fact.Y (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
367 break; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
368 |
20762
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
369 case 2: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
370 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
371 PermMatrix P = fact.P (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
372 FloatMatrix L = P.transpose () * fact.L (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
373 retval(1) = get_lu_u (fact); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
374 retval(0) = L; |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
375 } |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
376 break; |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
377 |
20762
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
378 case 3: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
379 default: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
380 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
381 if (vecout) |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
382 retval(2) = fact.P_vec (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
383 else |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
384 retval(2) = fact.P (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
385 retval(1) = get_lu_u (fact); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
386 retval(0) = get_lu_l (fact); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
387 } |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
388 break; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
389 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
390 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
391 else |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
392 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
393 Matrix m = arg.matrix_value (); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
394 |
20762
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
395 LU fact (m); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
396 |
20762
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
397 switch (nargout) |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
398 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
399 case 0: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
400 case 1: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
401 retval(0) = fact.Y (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
402 break; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
403 |
20762
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
404 case 2: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
405 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
406 PermMatrix P = fact.P (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
407 Matrix L = P.transpose () * fact.L (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
408 retval(1) = get_lu_u (fact); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
409 retval(0) = L; |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
410 } |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
411 break; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
412 |
20762
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
413 case 3: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
414 default: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
415 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
416 if (vecout) |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
417 retval(2) = fact.P_vec (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
418 else |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
419 retval(2) = fact.P (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
420 retval(1) = get_lu_u (fact); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
421 retval(0) = get_lu_l (fact); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
422 } |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
423 break; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
424 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
425 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
426 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
427 else if (arg.is_complex_type ()) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
428 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
429 if (arg.is_single_type ()) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
430 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
431 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
|
432 |
20762
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
433 FloatComplexLU fact (m); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
434 |
20762
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
435 switch (nargout) |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
436 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
437 case 0: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
438 case 1: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
439 retval(0) = fact.Y (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
440 break; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
441 |
20762
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
442 case 2: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
443 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
444 PermMatrix P = fact.P (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
445 FloatComplexMatrix L = P.transpose () * fact.L (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
446 retval(1) = get_lu_u (fact); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
447 retval(0) = L; |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
448 } |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
449 break; |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
450 |
20762
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
451 case 3: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
452 default: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
453 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
454 if (vecout) |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
455 retval(2) = fact.P_vec (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
456 else |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
457 retval(2) = fact.P (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
458 retval(1) = get_lu_u (fact); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
459 retval(0) = get_lu_l (fact); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
460 } |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
461 break; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
462 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
463 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
464 else |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
465 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
466 ComplexMatrix m = arg.complex_matrix_value (); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
467 |
20762
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
468 ComplexLU fact (m); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
469 |
20762
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
470 switch (nargout) |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
471 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
472 case 0: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
473 case 1: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
474 retval(0) = fact.Y (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
475 break; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
476 |
20762
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
477 case 2: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
478 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
479 PermMatrix P = fact.P (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
480 ComplexMatrix L = P.transpose () * fact.L (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
481 retval(1) = get_lu_u (fact); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
482 retval(0) = L; |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
483 } |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
484 break; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
485 |
20762
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
486 case 3: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
487 default: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
488 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
489 if (vecout) |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
490 retval(2) = fact.P_vec (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
491 else |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
492 retval(2) = fact.P (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
493 retval(1) = get_lu_u (fact); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
494 retval(0) = get_lu_l (fact); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
495 } |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20438
diff
changeset
|
496 break; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
497 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
498 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
499 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
500 else |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
501 gripe_wrong_type_arg ("lu", arg); |
2928 | 502 } |
503 | |
504 return retval; | |
505 } | |
506 | |
507 /* | |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
508 %!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
|
509 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
510 %!test |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
511 %! [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
|
512 %! 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
|
513 %! 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
|
514 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
515 %!test |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
516 %! [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
|
517 %! 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
|
518 %! 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
|
519 %! 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
|
520 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
521 %!test |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
522 %! [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
|
523 %! 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
|
524 %! 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
|
525 %! 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
|
526 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
527 %!test |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
528 %! [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
|
529 %! 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
|
530 %! 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
|
531 %! 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
|
532 |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
533 %!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
|
534 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
535 %!test |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
536 %! [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
|
537 %! 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
|
538 %! 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
|
539 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
540 %!test |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
541 %! [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
|
542 %! 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
|
543 %! 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
|
544 %! 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
|
545 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
546 %!test |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
547 %! [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
|
548 %! 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
|
549 %! 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
|
550 %! 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
|
551 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
552 %!test |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
553 %! [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
|
554 %! 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
|
555 %! 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
|
556 %! 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
|
557 |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
558 %!error lu () |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
559 %!error <can not define pivoting threshold> lu ([1, 2; 3, 4], 2) |
18493
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
560 |
18556
4c2465444a96
Use %!testif HAVE_UMFPACK for sparse lu test added in cset 2a45b6b87bee
Mike Miller <mtmiller@ieee.org>
parents:
18495
diff
changeset
|
561 %!testif HAVE_UMFPACK |
18493
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
562 %! 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
|
563 %! 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
|
564 %! Bv = [1 1 1 1 1 1 -1 1 1 1 1 -1 1 -1 1]; |
18495
93c019f00e59
maint: whitespace fixes for libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
18494
diff
changeset
|
565 %! B = sparse (Bi, Bj, Bv); |
93c019f00e59
maint: whitespace fixes for libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
18494
diff
changeset
|
566 %! [L, U] = lu (B); |
18494
807e21ebddd5
maint: Style fixes for tests in libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
18493
diff
changeset
|
567 %! assert (L*U, B); |
18495
93c019f00e59
maint: whitespace fixes for libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
18494
diff
changeset
|
568 %! [L, U, P] = lu(B); |
18494
807e21ebddd5
maint: Style fixes for tests in libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
18493
diff
changeset
|
569 %! assert (P'*L*U, B); |
18495
93c019f00e59
maint: whitespace fixes for libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
18494
diff
changeset
|
570 %! [L, U, P, Q] = lu (B); |
18494
807e21ebddd5
maint: Style fixes for tests in libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
18493
diff
changeset
|
571 %! assert (P'*L*U*Q', B); |
18493
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
572 |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
573 */ |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
574 |
9708 | 575 static |
576 bool check_lu_dims (const octave_value& l, const octave_value& u, | |
577 const octave_value& p) | |
578 { | |
18099
6a71e5030df5
Follow coding convention of defining and initializing only 1 variable per line in liboctinterp.
Rik <rik@octave.org>
parents:
17787
diff
changeset
|
579 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
|
580 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
|
581 octave_idx_type n = u.columns (); |
9708 | 582 return ((l.ndims () == 2 && u.ndims () == 2 && k == l.columns ()) |
20071
17d647821d61
maint: More cleanup of C++ code to follow Octave coding conventions.
John W. Eaton <jwe@octave.org>
parents:
20068
diff
changeset
|
583 && k == std::min (m, n) |
17d647821d61
maint: More cleanup of C++ code to follow Octave coding conventions.
John W. Eaton <jwe@octave.org>
parents:
20068
diff
changeset
|
584 && (p.is_undefined () || p.rows () == m)); |
9708 | 585 } |
586 | |
15039
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14854
diff
changeset
|
587 DEFUN (luupdate, args, , |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
588 "-*- texinfo -*-\n\ |
15039
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14854
diff
changeset
|
589 @deftypefn {Built-in Function} {[@var{L}, @var{U}] =} luupdate (@var{L}, @var{U}, @var{x}, @var{y})\n\ |
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14854
diff
changeset
|
590 @deftypefnx {Built-in Function} {[@var{L}, @var{U}, @var{P}] =} luupdate (@var{L}, @var{U}, @var{P}, @var{x}, @var{y})\n\ |
9709 | 591 Given an LU@tie{}factorization of a real or complex matrix\n\ |
592 @w{@var{A} = @var{L}*@var{U}}, @var{L}@tie{}lower unit trapezoidal and\n\ | |
593 @var{U}@tie{}upper trapezoidal, return the LU@tie{}factorization\n\ | |
594 of @w{@var{A} + @var{x}*@var{y}.'}, where @var{x} and @var{y} are\n\ | |
595 column vectors (rank-1 update) or matrices with equal number of columns\n\ | |
596 (rank-k update).\n\ | |
20382
4f45eaf83908
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
20071
diff
changeset
|
597 \n\ |
4f45eaf83908
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
20071
diff
changeset
|
598 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:
20071
diff
changeset
|
599 (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:
20071
diff
changeset
|
600 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:
20071
diff
changeset
|
601 LU@tie{}factorization as obtained by @code{lu}:\n\ |
9709 | 602 \n\ |
603 @example\n\ | |
14360
97883071e8e4
doc: Correct off-by-1 spacings in all .cc docstrings
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
604 [@var{L}, @var{U}, @var{P}] = lu (@var{A});\n\ |
9709 | 605 @end example\n\ |
606 \n\ | |
10846
a4f482e66b65
Grammarcheck more of the documentation.
Rik <octave@nomad.inbox5.com>
parents:
10840
diff
changeset
|
607 @noindent\n\ |
17268
1c21f264d26f
doc: Rename @xcode macro to @tcode (transpose code) for clarity.
Rik <rik@octave.org>
parents:
16920
diff
changeset
|
608 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
|
609 either as\n\ |
9709 | 610 \n\ |
611 @example\n\ | |
14360
97883071e8e4
doc: Correct off-by-1 spacings in all .cc docstrings
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
612 [@var{L1}, @var{U1}] = lu (@var{L}, @var{U}, @var{P}*@var{x}, @var{y})\n\ |
9709 | 613 @end example\n\ |
614 \n\ | |
10846
a4f482e66b65
Grammarcheck more of the documentation.
Rik <octave@nomad.inbox5.com>
parents:
10840
diff
changeset
|
615 @noindent\n\ |
9709 | 616 or\n\ |
617 \n\ | |
618 @example\n\ | |
14360
97883071e8e4
doc: Correct off-by-1 spacings in all .cc docstrings
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
619 [@var{L1}, @var{U1}, @var{P1}] = lu (@var{L}, @var{U}, @var{P}, @var{x}, @var{y})\n\ |
9709 | 620 @end example\n\ |
621 \n\ | |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
622 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
|
623 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
|
624 stable.\n\ |
9709 | 625 \n\ |
20382
4f45eaf83908
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
20071
diff
changeset
|
626 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:
20071
diff
changeset
|
627 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:
20071
diff
changeset
|
628 factorization from scratch.\n\ |
16920
53eaa83e4181
doc: Add seealso links between various factorization forms.
Rik <rik@octave.org>
parents:
15195
diff
changeset
|
629 @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
|
630 @end deftypefn") |
9708 | 631 { |
632 octave_idx_type nargin = args.length (); | |
633 octave_value_list retval; | |
634 | |
635 bool pivoted = nargin == 5; | |
636 | |
637 if (nargin != 4 && nargin != 5) | |
638 { | |
639 print_usage (); | |
640 return retval; | |
641 } | |
642 | |
643 octave_value argl = args(0); | |
644 octave_value argu = args(1); | |
645 octave_value argp = pivoted ? args(2) : octave_value (); | |
646 octave_value argx = args(2 + pivoted); | |
647 octave_value argy = args(3 + pivoted); | |
648 | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
649 if (argl.is_numeric_type () && argu.is_numeric_type () |
9708 | 650 && argx.is_numeric_type () && argy.is_numeric_type () |
651 && (! pivoted || argp.is_perm_matrix ())) | |
652 { | |
653 if (check_lu_dims (argl, argu, argp)) | |
654 { | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
655 PermMatrix P = (pivoted |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
656 ? argp.perm_matrix_value () |
9708 | 657 : PermMatrix::eye (argl.rows ())); |
658 | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
659 if (argl.is_real_type () |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
660 && argu.is_real_type () |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
661 && argx.is_real_type () |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
662 && argy.is_real_type ()) |
9708 | 663 { |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
664 // all real case |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
665 if (argl.is_single_type () |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
666 || argu.is_single_type () |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
667 || argx.is_single_type () |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
668 || argy.is_single_type ()) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
669 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
670 FloatMatrix L = argl.float_matrix_value (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
671 FloatMatrix U = argu.float_matrix_value (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
672 FloatMatrix x = argx.float_matrix_value (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
673 FloatMatrix y = argy.float_matrix_value (); |
9708 | 674 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
675 FloatLU fact (L, U, P); |
9708 | 676 if (pivoted) |
677 fact.update_piv (x, y); | |
678 else | |
679 fact.update (x, y); | |
680 | |
681 if (pivoted) | |
682 retval(2) = fact.P (); | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
683 retval(1) = get_lu_u (fact); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
684 retval(0) = get_lu_l (fact); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
685 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
686 else |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
687 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
688 Matrix L = argl.matrix_value (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
689 Matrix U = argu.matrix_value (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
690 Matrix x = argx.matrix_value (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
691 Matrix y = argy.matrix_value (); |
9708 | 692 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
693 LU fact (L, U, P); |
9708 | 694 if (pivoted) |
695 fact.update_piv (x, y); | |
696 else | |
697 fact.update (x, y); | |
698 | |
699 if (pivoted) | |
700 retval(2) = fact.P (); | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
701 retval(1) = get_lu_u (fact); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
702 retval(0) = get_lu_l (fact); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
703 } |
9708 | 704 } |
705 else | |
706 { | |
707 // complex case | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
708 if (argl.is_single_type () |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
709 || argu.is_single_type () |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
710 || argx.is_single_type () |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
711 || argy.is_single_type ()) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
712 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
713 FloatComplexMatrix L = argl.float_complex_matrix_value (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
714 FloatComplexMatrix U = argu.float_complex_matrix_value (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
715 FloatComplexMatrix x = argx.float_complex_matrix_value (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
716 FloatComplexMatrix y = argy.float_complex_matrix_value (); |
9708 | 717 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
718 FloatComplexLU fact (L, U, P); |
9708 | 719 if (pivoted) |
720 fact.update_piv (x, y); | |
721 else | |
722 fact.update (x, y); | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
723 |
9708 | 724 if (pivoted) |
725 retval(2) = fact.P (); | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
726 retval(1) = get_lu_u (fact); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
727 retval(0) = get_lu_l (fact); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
728 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
729 else |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
730 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
731 ComplexMatrix L = argl.complex_matrix_value (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
732 ComplexMatrix U = argu.complex_matrix_value (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
733 ComplexMatrix x = argx.complex_matrix_value (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
734 ComplexMatrix y = argy.complex_matrix_value (); |
9708 | 735 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
736 ComplexLU fact (L, U, P); |
9708 | 737 if (pivoted) |
738 fact.update_piv (x, y); | |
739 else | |
740 fact.update (x, y); | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
741 |
9708 | 742 if (pivoted) |
743 retval(2) = fact.P (); | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
744 retval(1) = get_lu_u (fact); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
745 retval(0) = get_lu_l (fact); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
746 } |
9708 | 747 } |
748 } | |
749 else | |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
750 error ("luupdate: dimension mismatch"); |
9708 | 751 } |
752 else | |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
753 error ("luupdate: L, U, X, and Y must be numeric"); |
9708 | 754 |
755 return retval; | |
756 } | |
757 | |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
758 /* |
10080
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
759 %!shared A, u, v, Ac, uc, vc |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
760 %! A = [0.091364 0.613038 0.999083; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
761 %! 0.594638 0.425302 0.603537; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
762 %! 0.383594 0.291238 0.085574; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
763 %! 0.265712 0.268003 0.238409; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
764 %! 0.669966 0.743851 0.445057 ]; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
765 %! |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
766 %! u = [0.85082; |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
767 %! 0.76426; |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
768 %! 0.42883; |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
769 %! 0.53010; |
10080
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
770 %! 0.80683 ]; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
771 %! |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
772 %! v = [0.98810; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
773 %! 0.24295; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
774 %! 0.43167 ]; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
775 %! |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
776 %! 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
|
777 %! 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
|
778 %! 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
|
779 %! 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
|
780 %! 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
|
781 %! |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
782 %! 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
|
783 %! 0.13141 + 0.43708i; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
784 %! 0.29808 + 0.08789i; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
785 %! 0.69821 + 0.38844i; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
786 %! 0.74871 + 0.25821i ]; |
10080
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
787 %! |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
788 %! 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
|
789 %! 0.20820 + 0.93090i; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
790 %! 0.86184 + 0.34689i ]; |
10080
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
791 %! |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
792 |
10094
ab1101011a6d
lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents:
10080
diff
changeset
|
793 %!testif HAVE_QRUPDATE_LUU |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
794 %! [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
|
795 %! [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
|
796 %! 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
|
797 %! 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
|
798 %! 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
|
799 %! |
10094
ab1101011a6d
lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents:
10080
diff
changeset
|
800 %!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
|
801 %! [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
|
802 %! [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
|
803 %! 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
|
804 %! 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
|
805 %! 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
|
806 |
10094
ab1101011a6d
lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents:
10080
diff
changeset
|
807 %!testif HAVE_QRUPDATE_LUU |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
808 %! [L,U,P] = lu (single (A)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
809 %! [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
|
810 %! assert (norm (vec (tril (L)-L), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
811 %! assert (norm (vec (triu (U)-U), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
812 %! assert (norm (vec (P'*L*U - 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
|
813 %! |
10094
ab1101011a6d
lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents:
10080
diff
changeset
|
814 %!testif HAVE_QRUPDATE_LUU |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
815 %! [L,U,P] = lu (single (Ac)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
816 %! [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
|
817 %! assert (norm (vec (tril (L)-L), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
818 %! assert (norm (vec (triu (U)-U), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
819 %! assert (norm (vec (P'*L*U - single (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
|
820 |
10094
ab1101011a6d
lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents:
10080
diff
changeset
|
821 %!testif HAVE_QRUPDATE_LUU |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
822 %! [L,U,P] = lu (A); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
823 %! [L,U,P] = luupdate (L,U,P,u,v); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
824 %! assert (norm (vec (tril (L)-L), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
825 %! assert (norm (vec (triu (U)-U), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
826 %! assert (norm (vec (P'*L*U - 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
|
827 %! |
10094
ab1101011a6d
lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents:
10080
diff
changeset
|
828 %!testif HAVE_QRUPDATE_LUU |
19039
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18833
diff
changeset
|
829 %! [L,U,P] = lu (A); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18833
diff
changeset
|
830 %! [~,ordcols] = max (P,[],1); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18833
diff
changeset
|
831 %! [~,ordrows] = max (P,[],2); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18833
diff
changeset
|
832 %! P1 = eye (size(P))(:,ordcols); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18833
diff
changeset
|
833 %! P2 = eye (size(P))(ordrows,:); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18833
diff
changeset
|
834 %! assert(P1 == P); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18833
diff
changeset
|
835 %! assert(P2 == P); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18833
diff
changeset
|
836 %! [L,U,P] = luupdate (L,U,P,u,v); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18833
diff
changeset
|
837 %! [L,U,P1] = luupdate (L,U,P1,u,v); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18833
diff
changeset
|
838 %! [L,U,P2] = luupdate (L,U,P2,u,v); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18833
diff
changeset
|
839 %! assert(P1 == P); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18833
diff
changeset
|
840 %! assert(P2 == P); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18833
diff
changeset
|
841 %! |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18833
diff
changeset
|
842 %!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
|
843 %! [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
|
844 %! [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
|
845 %! 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
|
846 %! 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
|
847 %! 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
|
848 |
10094
ab1101011a6d
lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents:
10080
diff
changeset
|
849 %!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
|
850 %! [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
|
851 %! [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
|
852 %! 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
|
853 %! 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
|
854 %! 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
|
855 %! |
10094
ab1101011a6d
lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents:
10080
diff
changeset
|
856 %!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
|
857 %! [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
|
858 %! [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
|
859 %! 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
|
860 %! 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
|
861 %! 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
|
862 */ |