Mercurial > hg > octave-nkf
annotate src/DLD-FUNCTIONS/lu.cc @ 8869:c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
author | Jaroslav Hajek <highegg@gmail.com> |
---|---|
date | Wed, 25 Feb 2009 09:10:38 +0100 |
parents | 445d27d79f4e |
children | eb63fbe60fab |
rev | line source |
---|---|
2928 | 1 /* |
2 | |
7017 | 3 Copyright (C) 1996, 1997, 1999, 2000, 2003, 2005, 2006, 2007 |
4 John W. Eaton | |
2928 | 5 |
6 This file is part of Octave. | |
7 | |
8 Octave is free software; you can redistribute it and/or modify it | |
9 under the terms of the GNU General Public License as published by the | |
7016 | 10 Free Software Foundation; either version 3 of the License, or (at your |
11 option) any later version. | |
2928 | 12 |
13 Octave is distributed in the hope that it will be useful, but WITHOUT | |
14 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or | |
15 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License | |
16 for more details. | |
17 | |
18 You should have received a copy of the GNU General Public License | |
7016 | 19 along with Octave; see the file COPYING. If not, see |
20 <http://www.gnu.org/licenses/>. | |
2928 | 21 |
22 */ | |
23 | |
24 #ifdef HAVE_CONFIG_H | |
25 #include <config.h> | |
26 #endif | |
27 | |
28 #include "CmplxLU.h" | |
29 #include "dbleLU.h" | |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
30 #include "fCmplxLU.h" |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
31 #include "floatLU.h" |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
32 #include "SparseCmplxLU.h" |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
33 #include "SparsedbleLU.h" |
2928 | 34 |
35 #include "defun-dld.h" | |
36 #include "error.h" | |
37 #include "gripes.h" | |
38 #include "oct-obj.h" | |
39 #include "utils.h" | |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
40 #include "ov-re-sparse.h" |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
41 #include "ov-cx-sparse.h" |
2928 | 42 |
8869
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
43 template <class MT> |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
44 static octave_value |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
45 maybe_set_triangular (const MT& m, MatrixType::matrix_type t = MatrixType::Upper) |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
46 { |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
47 typedef typename MT::element_type T; |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
48 octave_value retval; |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
49 octave_idx_type r = m.rows (), c = m.columns (); |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
50 if (r == c) |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
51 { |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
52 const T zero = T(); |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
53 octave_idx_type i = 0; |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
54 for (;i != r && m(i,i) != zero; i++) ; |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
55 if (i == r) |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
56 retval = octave_value (m, MatrixType (t)); |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
57 else |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
58 retval = m; |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
59 } |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
60 else |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
61 retval = m; |
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 return retval; |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
64 } |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
65 |
2928 | 66 DEFUN_DLD (lu, args, nargout, |
3548 | 67 "-*- texinfo -*-\n\ |
3372 | 68 @deftypefn {Loadable Function} {[@var{l}, @var{u}, @var{p}] =} lu (@var{a})\n\ |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
69 @deftypefnx {Loadable Function} {[@var{l}, @var{u}, @var{p}, @var{q}] =} lu (@var{s})\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
70 @deftypefnx {Loadable Function} {[@var{l}, @var{u}, @var{p}, @var{q}, @var{r}] =} lu (@var{s})\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
71 @deftypefnx {Loadable Function} {[@dots{}] =} lu (@var{s}, @var{thres})\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
72 @deftypefnx {Loadable Function} {@var{y} =} lu (@dots{})\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
73 @deftypefnx {Loadable Function} {[@dots{}] =} lu (@dots{}, 'vector')\n\ |
3372 | 74 @cindex LU decomposition\n\ |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
75 Compute the LU decomposition of @var{a}. If @var{a} is full subroutines from\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
76 @sc{Lapack} are used and if @var{a} is sparse then UMFPACK is used. The\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
77 result is returned in a permuted form, according to the optional return\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
78 value @var{p}. For example, given the matrix @code{a = [1, 2; 3, 4]},\n\ |
3372 | 79 \n\ |
80 @example\n\ | |
81 [l, u, p] = lu (a)\n\ | |
82 @end example\n\ | |
83 \n\ | |
84 @noindent\n\ | |
85 returns\n\ | |
86 \n\ | |
87 @example\n\ | |
88 l =\n\ | |
89 \n\ | |
90 1.00000 0.00000\n\ | |
91 0.33333 1.00000\n\ | |
92 \n\ | |
93 u =\n\ | |
94 \n\ | |
95 3.00000 4.00000\n\ | |
96 0.00000 0.66667\n\ | |
97 \n\ | |
98 p =\n\ | |
99 \n\ | |
100 0 1\n\ | |
101 1 0\n\ | |
102 @end example\n\ | |
4329 | 103 \n\ |
5457 | 104 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
|
105 \n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
106 Called with two or three output arguments and a spare input matrix,\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
107 then @dfn{lu} does not attempt to perform sparsity preserving column\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
108 permutations. Called with a fourth output argument, the sparsity\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
109 preserving column transformation @var{Q} is returned, such that\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
110 @code{@var{p} * @var{a} * @var{q} = @var{l} * @var{u}}.\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
111 \n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
112 Called with a fifth output argument and a sparse input matrix, then\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
113 @dfn{lu} attempts to use a scaling factor @var{r} on the input matrix\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
114 such that @code{@var{p} * (@var{r} \\ @var{a}) * @var{q} = @var{l} * @var{u}}.\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
115 This typically leads to a sparser and more stable factorsation.\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
116 \n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
117 An additional input argument @var{thres}, that defines the pivoting\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
118 threshold can be given. @var{thres} can be a scalar, in which case\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
119 it defines UMFPACK pivoting tolerance for both symmetric and unsymmetric\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
120 cases. If @var{thres} is a two element vector, then the first element\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
121 defines the pivoting tolerance for the unsymmetric UMFPACK pivoting\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
122 strategy and the second the symmetric strategy. By default, the values\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
123 defined by @code{spparms} are used and are by default @code{[0.1, 0.001]}.\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
124 \n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
125 Given the string argument 'vector', @dfn{lu} returns the values of @var{p}\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
126 @var{q} as vector values, such that for full matrix, @code{@var{a}\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
127 (@var{p},:) = @var{l} * @var{u}}, and @code{@var{r}(@var{p},:) * @var{a}\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
128 (:, @var{q}) = @var{l} * @var{u}}.\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
129 \n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
130 With two output arguments, returns the permuted forms of the upper and\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
131 lower triangular matrices, such that @code{@var{a} = @var{l} * @var{u}}.\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
132 With one output argument @var{y}, then the matrix returned by the @sc{Lapack}\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
133 routines is returned. If the input matrix is sparse then the matrix @var{l}\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
134 is embedded into @var{u} to give a return value similar to the full case.\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
135 For both full and sparse matrices, @dfn{lu} looses the permutation\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
136 information.\n\ |
3372 | 137 @end deftypefn") |
2928 | 138 { |
139 octave_value_list retval; | |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
140 int nargin = args.length (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
141 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
|
142 bool scale = (nargout == 5); |
2928 | 143 |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
144 if (nargin < 1 || (issparse && (nargin > 3 || nargout > 5)) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
145 || (!issparse && (nargin > 2 || nargout > 3))) |
2928 | 146 { |
5823 | 147 print_usage (); |
2928 | 148 return retval; |
149 } | |
150 | |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
151 bool vecout = false; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
152 Matrix thres; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
153 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
154 int n = 1; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
155 while (n < nargin && ! error_state) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
156 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
157 if (args (n).is_string ()) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
158 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
159 std::string tmp = args(n++).string_value (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
160 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
161 if (! error_state ) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
162 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
163 if (tmp.compare ("vector") == 0) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
164 vecout = true; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
165 else |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
166 error ("lu: unrecognized string argument"); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
167 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
168 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
169 else |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
170 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
171 Matrix tmp = args(n++).matrix_value (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
172 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
173 if (! error_state ) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
174 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
175 if (!issparse) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
176 error ("lu: can not define pivoting threshold for full matrices"); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
177 else if (tmp.nelem () == 1) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
178 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
179 thres.resize(1,2); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
180 thres(0) = tmp(0); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
181 thres(1) = tmp(0); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
182 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
183 else if (tmp.nelem () == 2) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
184 thres = tmp; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
185 else |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
186 error ("lu: expecting 2 element vector for thres"); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
187 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
188 } |
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) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
201 return retval; |
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) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
203 return octave_value_list (5, SparseMatrix ()); |
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 ColumnVector Qinit; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
206 if (nargout < 4) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
207 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
208 Qinit.resize (nc); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
209 for (octave_idx_type i = 0; i < nc; i++) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
210 Qinit (i) = i; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
211 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
212 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
213 if (arg.is_real_type ()) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
214 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
215 SparseMatrix m = arg.sparse_matrix_value (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
216 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
217 switch (nargout) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
218 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
219 case 0: |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
220 case 1: |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
221 case 2: |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
222 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
223 SparseLU fact (m, Qinit, thres, false, true); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
224 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
225 if (nargout < 2) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
226 retval (0) = fact.Y (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
227 else |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
228 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
229 SparseMatrix P = fact.Pr (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
230 SparseMatrix L = P.transpose () * fact.L (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
231 retval(1) = octave_value (fact.U (), |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
232 MatrixType (MatrixType::Upper)); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
233 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
234 retval(0) = octave_value (L, |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
235 MatrixType (MatrixType::Permuted_Lower, |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
236 nr, fact.row_perm ())); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
237 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
238 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
239 break; |
2928 | 240 |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
241 case 3: |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
242 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
243 SparseLU fact (m, Qinit, thres, false, true); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
244 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
245 if (vecout) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
246 retval (2) = fact.Pr_vec (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
247 else |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
248 retval(2) = fact.Pr (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
249 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
250 retval(1) = octave_value (fact.U (), |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
251 MatrixType (MatrixType::Upper)); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
252 retval(0) = octave_value (fact.L (), |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
253 MatrixType (MatrixType::Lower)); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
254 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
255 break; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
256 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
257 case 4: |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
258 default: |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
259 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
260 SparseLU fact (m, thres, scale); |
2928 | 261 |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
262 if (scale) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
263 retval(4) = fact.R (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
264 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
265 if (vecout) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
266 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
267 retval(3) = fact.Pc_vec (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
268 retval(2) = fact.Pr_vec (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
269 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
270 else |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
271 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
272 retval(3) = fact.Pc (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
273 retval(2) = fact.Pr (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
274 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
275 retval(1) = octave_value (fact.U (), |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
276 MatrixType (MatrixType::Upper)); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
277 retval(0) = octave_value (fact.L (), |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
278 MatrixType (MatrixType::Lower)); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
279 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
280 break; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
281 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
282 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
283 else if (arg.is_complex_type ()) |
2928 | 284 { |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
285 SparseComplexMatrix m = arg.sparse_complex_matrix_value (); |
2928 | 286 |
287 switch (nargout) | |
288 { | |
289 case 0: | |
290 case 1: | |
291 case 2: | |
292 { | |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
293 SparseComplexLU fact (m, Qinit, thres, false, true); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
294 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
295 if (nargout < 2) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
296 retval (0) = fact.Y (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
297 else |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
298 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
299 SparseMatrix P = fact.Pr (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
300 SparseComplexMatrix L = P.transpose () * fact.L (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
301 retval(1) = octave_value (fact.U (), |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
302 MatrixType (MatrixType::Upper)); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
303 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
304 retval(0) = octave_value (L, |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
305 MatrixType (MatrixType::Permuted_Lower, |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
306 nr, fact.row_perm ())); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
307 } |
2928 | 308 } |
309 break; | |
310 | |
311 case 3: | |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
312 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
313 SparseComplexLU fact (m, Qinit, thres, false, true); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
314 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
315 if (vecout) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
316 retval (2) = fact.Pr_vec (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
317 else |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
318 retval(2) = fact.Pr (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
319 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
320 retval(1) = octave_value (fact.U (), |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
321 MatrixType (MatrixType::Upper)); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
322 retval(0) = octave_value (fact.L (), |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
323 MatrixType (MatrixType::Lower)); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
324 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
325 break; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
326 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
327 case 4: |
2928 | 328 default: |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
329 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
330 SparseComplexLU fact (m, thres, scale); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
331 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
332 if (scale) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
333 retval(4) = fact.R (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
334 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
335 if (vecout) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
336 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
337 retval(3) = fact.Pc_vec (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
338 retval(2) = fact.Pr_vec (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
339 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
340 else |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
341 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
342 retval(3) = fact.Pc (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
343 retval(2) = fact.Pr (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
344 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
345 retval(1) = octave_value (fact.U (), |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
346 MatrixType (MatrixType::Upper)); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
347 retval(0) = octave_value (fact.L (), |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
348 MatrixType (MatrixType::Lower)); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
349 } |
2928 | 350 break; |
351 } | |
352 } | |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
353 else |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
354 gripe_wrong_type_arg ("lu", arg); |
2928 | 355 } |
356 else | |
357 { | |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
358 if (arg_is_empty < 0) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
359 return retval; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
360 else if (arg_is_empty > 0) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
361 return octave_value_list (3, Matrix ()); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
362 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
363 if (arg.is_real_type ()) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
364 { |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
365 if (arg.is_single_type ()) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
366 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
367 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
|
368 |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
369 if (! error_state) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
370 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
371 FloatLU fact (m); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
372 |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
373 switch (nargout) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
374 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
375 case 0: |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
376 case 1: |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
377 retval(0) = fact.Y (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
378 break; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
379 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
380 case 2: |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
381 { |
8367
445d27d79f4e
support permutation matrix objects
Jaroslav Hajek <highegg@gmail.com>
parents:
7814
diff
changeset
|
382 PermMatrix P = fact.P (); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
383 FloatMatrix L = P.transpose () * fact.L (); |
8869
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
384 retval(1) = maybe_set_triangular (fact.U (), MatrixType::Upper); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
385 retval(0) = L; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
386 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
387 break; |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
388 |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
389 case 3: |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
390 default: |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
391 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
392 if (vecout) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
393 retval(2) = fact.P_vec (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
394 else |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
395 retval(2) = fact.P (); |
8869
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
396 retval(1) = maybe_set_triangular (fact.U (), MatrixType::Upper); |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
397 retval(0) = maybe_set_triangular (fact.L (), MatrixType::Lower); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
398 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
399 break; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
400 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
401 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
402 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
403 else |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
404 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
405 Matrix m = arg.matrix_value (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
406 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
407 if (! error_state) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
408 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
409 LU fact (m); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
410 |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
411 switch (nargout) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
412 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
413 case 0: |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
414 case 1: |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
415 retval(0) = fact.Y (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
416 break; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
417 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
418 case 2: |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
419 { |
8367
445d27d79f4e
support permutation matrix objects
Jaroslav Hajek <highegg@gmail.com>
parents:
7814
diff
changeset
|
420 PermMatrix P = fact.P (); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
421 Matrix L = P.transpose () * fact.L (); |
8869
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
422 retval(1) = maybe_set_triangular (fact.U (), MatrixType::Upper); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
423 retval(0) = L; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
424 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
425 break; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
426 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
427 case 3: |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
428 default: |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
429 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
430 if (vecout) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
431 retval(2) = fact.P_vec (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
432 else |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
433 retval(2) = fact.P (); |
8869
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
434 retval(1) = maybe_set_triangular (fact.U (), MatrixType::Upper); |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
435 retval(0) = maybe_set_triangular (fact.L (), MatrixType::Lower); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
436 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
437 break; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
438 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
439 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
440 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
441 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
442 else if (arg.is_complex_type ()) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
443 { |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
444 if (arg.is_single_type ()) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
445 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
446 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
|
447 |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
448 if (! error_state) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
449 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
450 FloatComplexLU fact (m); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
451 |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
452 switch (nargout) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
453 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
454 case 0: |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
455 case 1: |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
456 retval(0) = fact.Y (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
457 break; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
458 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
459 case 2: |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
460 { |
8367
445d27d79f4e
support permutation matrix objects
Jaroslav Hajek <highegg@gmail.com>
parents:
7814
diff
changeset
|
461 PermMatrix P = fact.P (); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
462 FloatComplexMatrix L = P.transpose () * fact.L (); |
8869
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
463 retval(1) = maybe_set_triangular (fact.U (), MatrixType::Upper); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
464 retval(0) = L; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
465 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
466 break; |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
467 |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
468 case 3: |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
469 default: |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
470 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
471 if (vecout) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
472 retval(2) = fact.P_vec (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
473 else |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
474 retval(2) = fact.P (); |
8869
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
475 retval(1) = maybe_set_triangular (fact.U (), MatrixType::Upper); |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
476 retval(0) = maybe_set_triangular (fact.L (), MatrixType::Lower); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
477 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
478 break; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
479 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
480 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
481 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
482 else |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
483 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
484 ComplexMatrix m = arg.complex_matrix_value (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
485 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
486 if (! error_state) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
487 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
488 ComplexLU fact (m); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
489 |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
490 switch (nargout) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
491 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
492 case 0: |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
493 case 1: |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
494 retval(0) = fact.Y (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
495 break; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
496 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
497 case 2: |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
498 { |
8367
445d27d79f4e
support permutation matrix objects
Jaroslav Hajek <highegg@gmail.com>
parents:
7814
diff
changeset
|
499 PermMatrix P = fact.P (); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
500 ComplexMatrix L = P.transpose () * fact.L (); |
8869
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
501 retval(1) = maybe_set_triangular (fact.U (), MatrixType::Upper); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
502 retval(0) = L; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
503 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
504 break; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
505 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
506 case 3: |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
507 default: |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
508 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
509 if (vecout) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
510 retval(2) = fact.P_vec (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
511 else |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
512 retval(2) = fact.P (); |
8869
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
513 retval(1) = maybe_set_triangular (fact.U (), MatrixType::Upper); |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
514 retval(0) = maybe_set_triangular (fact.L (), MatrixType::Lower); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
515 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
516 break; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
517 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
518 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
519 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
520 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
521 else |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
522 gripe_wrong_type_arg ("lu", arg); |
2928 | 523 } |
524 | |
525 return retval; | |
526 } | |
527 | |
528 /* | |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
529 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
530 %!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
|
531 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
532 %!test |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
533 %! [l, u] = lu ([1, 2; 3, 4]); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
534 %! assert(l, [1/3, 1; 1, 0], sqrt (eps)); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
535 %! assert(u, [3, 4; 0, 2/3], sqrt (eps)); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
536 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
537 %!test |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
538 %! [l, u, p] = lu ([1, 2; 3, 4]); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
539 %! assert(l, [1, 0; 1/3, 1], sqrt (eps)); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
540 %! assert(u, [3, 4; 0, 2/3], sqrt (eps)); |
8367
445d27d79f4e
support permutation matrix objects
Jaroslav Hajek <highegg@gmail.com>
parents:
7814
diff
changeset
|
541 %! 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
|
542 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
543 %!test |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
544 %! [l, u, p] = lu ([1, 2; 3, 4],'vector'); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
545 %! assert(l, [1, 0; 1/3, 1], sqrt (eps)); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
546 %! assert(u, [3, 4; 0, 2/3], sqrt (eps)); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
547 %! assert(p, [2;1], sqrt (eps)); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
548 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
549 %!test |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
550 %! [l u p] = lu ([1, 2; 3, 4; 5, 6]); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
551 %! assert(l, [1, 0; 1/5, 1; 3/5, 1/2], sqrt (eps)); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
552 %! assert(u, [5, 6; 0, 4/5], sqrt (eps)); |
8367
445d27d79f4e
support permutation matrix objects
Jaroslav Hajek <highegg@gmail.com>
parents:
7814
diff
changeset
|
553 %! 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
|
554 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
555 %!assert(lu (single([1, 2; 3, 4])), single([3, 4; 1/3, 2/3]), eps('single')); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
556 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
557 %!test |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
558 %! [l, u] = lu (single([1, 2; 3, 4])); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
559 %! assert(l, single([1/3, 1; 1, 0]), sqrt (eps('single'))); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
560 %! assert(u, single([3, 4; 0, 2/3]), sqrt (eps('single'))); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
561 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
562 %!test |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
563 %! [l, u, p] = lu (single([1, 2; 3, 4])); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
564 %! assert(l, single([1, 0; 1/3, 1]), sqrt (eps('single'))); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
565 %! assert(u, single([3, 4; 0, 2/3]), sqrt (eps('single'))); |
8367
445d27d79f4e
support permutation matrix objects
Jaroslav Hajek <highegg@gmail.com>
parents:
7814
diff
changeset
|
566 %! 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
|
567 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
568 %!test |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
569 %! [l, u, p] = lu (single([1, 2; 3, 4]),'vector'); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
570 %! assert(l, single([1, 0; 1/3, 1]), sqrt (eps('single'))); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
571 %! assert(u, single([3, 4; 0, 2/3]), sqrt (eps('single'))); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
572 %! assert(p, single([2;1]), sqrt (eps('single'))); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
573 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
574 %!test |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
575 %! [l u p] = lu (single([1, 2; 3, 4; 5, 6])); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
576 %! assert(l, single([1, 0; 1/5, 1; 3/5, 1/2]), sqrt (eps('single'))); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
577 %! assert(u, single([5, 6; 0, 4/5]), sqrt (eps('single'))); |
8367
445d27d79f4e
support permutation matrix objects
Jaroslav Hajek <highegg@gmail.com>
parents:
7814
diff
changeset
|
578 %! 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
|
579 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
580 %!error <Invalid call to lu.*> lu (); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
581 %!error lu ([1, 2; 3, 4], 2); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
582 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
583 */ |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
584 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
585 /* |
2928 | 586 ;;; Local Variables: *** |
587 ;;; mode: C++ *** | |
588 ;;; End: *** | |
589 */ |