Mercurial > octave-nkf
annotate libinterp/corefcn/luinc.cc @ 20654: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 | 4bed806ee3d4 |
children |
rev | line source |
---|---|
5282 | 1 /* |
2 | |
19731
4197fc428c7d
maint: Update copyright notices for 2015.
John W. Eaton <jwe@octave.org>
parents:
19182
diff
changeset
|
3 Copyright (C) 2005-2015 David Bateman |
5282 | 4 |
7016 | 5 This file is part of Octave. |
6 | |
5282 | 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. | |
5282 | 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/>. | |
5282 | 20 |
21 */ | |
22 | |
23 #ifdef HAVE_CONFIG_H | |
24 #include <config.h> | |
25 #endif | |
26 | |
15039
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14854
diff
changeset
|
27 #include "defun.h" |
5282 | 28 #include "error.h" |
29 #include "gripes.h" | |
30 #include "oct-obj.h" | |
31 #include "utils.h" | |
32 #include "oct-map.h" | |
33 | |
5785 | 34 #include "MatrixType.h" |
5282 | 35 #include "SparseCmplxLU.h" |
36 #include "SparsedbleLU.h" | |
37 #include "ov-re-sparse.h" | |
38 #include "ov-cx-sparse.h" | |
39 | |
20207
4f45eaf83908
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19731
diff
changeset
|
40 // FIXME: Deprecated in 4.0 and should be removed in 4.4. |
19182 | 41 DEFUN (__luinc__, args, nargout, |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
42 "-*- texinfo -*-\n\ |
15039
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14854
diff
changeset
|
43 @deftypefn {Built-in Function} {[@var{L}, @var{U}, @var{P}, @var{Q}] =} luinc (@var{A}, '0')\n\ |
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14854
diff
changeset
|
44 @deftypefnx {Built-in Function} {[@var{L}, @var{U}, @var{P}, @var{Q}] =} luinc (@var{A}, @var{droptol})\n\ |
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14854
diff
changeset
|
45 @deftypefnx {Built-in Function} {[@var{L}, @var{U}, @var{P}, @var{Q}] =} luinc (@var{A}, @var{opts})\n\ |
5282 | 46 @cindex LU decomposition\n\ |
11593
1577c6f80926
Use non-breaking spaces between certain adjectives and their nouns in docstrings.
Rik <octave@nomad.inbox5.com>
parents:
11586
diff
changeset
|
47 Produce the incomplete LU@tie{}factorization of the sparse matrix @var{A}.\n\ |
20207
4f45eaf83908
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19731
diff
changeset
|
48 \n\ |
5282 | 49 Two types of incomplete factorization are possible, and the type\n\ |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
50 is determined by the second argument to @code{luinc}.\n\ |
5282 | 51 \n\ |
17281
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
15195
diff
changeset
|
52 Called with a second argument of @qcode{'0'}, the zero-level incomplete\n\ |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
53 LU@tie{}factorization is produced. This creates a factorization of @var{A}\n\ |
18851
9ac2357f19bc
doc: Replace "non-zero" with "nonzero" to match existing usage.
Rik <rik@octave.org>
parents:
18712
diff
changeset
|
54 where the position of the nonzero arguments correspond to the same\n\ |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
55 positions as in the matrix @var{A}.\n\ |
5282 | 56 \n\ |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
57 Alternatively, the fill-in of the incomplete LU@tie{}factorization can\n\ |
5282 | 58 be controlled through the variable @var{droptol} or the structure\n\ |
10711
fbd7843974fa
Periodic grammar check of documentation files to ensure common format.
Rik <octave@nomad.inbox5.com>
parents:
10155
diff
changeset
|
59 @var{opts}. The @sc{umfpack} multifrontal factorization code by Tim A.\n\ |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
60 Davis is used for the incomplete LU@tie{}factorization, (availability\n\ |
5282 | 61 @url{http://www.cise.ufl.edu/research/sparse/umfpack/})\n\ |
62 \n\ | |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
63 @var{droptol} determines the values below which the values in the\n\ |
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
64 LU@tie{} factorization are dropped and replaced by zero. It must be a\n\ |
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
65 positive scalar, and any values in the factorization whose absolute value\n\ |
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
66 are less than this value are dropped, expect if leaving them increase the\n\ |
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
67 sparsity of the matrix. Setting @var{droptol} to zero results in a complete\n\ |
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
68 LU@tie{}factorization which is the default.\n\ |
5282 | 69 \n\ |
70 @var{opts} is a structure containing one or more of the fields\n\ | |
71 \n\ | |
72 @table @code\n\ | |
73 @item droptol\n\ | |
9064
7c02ec148a3c
Check grammar on all .cc files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
74 The drop tolerance as above. If @var{opts} only contains @code{droptol}\n\ |
5282 | 75 then this is equivalent to using the variable @var{droptol}.\n\ |
76 \n\ | |
77 @item milu\n\ | |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
78 A logical variable flagging whether to use the modified incomplete\n\ |
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
79 LU@tie{} factorization. In the case that @code{milu} is true, the dropped\n\ |
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
80 values are subtracted from the diagonal of the matrix @var{U} of the\n\ |
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
81 factorization. The default is @code{false}.\n\ |
5282 | 82 \n\ |
83 @item udiag\n\ | |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
84 A logical variable that flags whether zero elements on the diagonal of\n\ |
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
85 @var{U} should be replaced with @var{droptol} to attempt to avoid singular\n\ |
9064
7c02ec148a3c
Check grammar on all .cc files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
86 factors. The default is @code{false}.\n\ |
5282 | 87 \n\ |
88 @item thresh\n\ | |
9064
7c02ec148a3c
Check grammar on all .cc files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
89 Defines the pivot threshold in the interval [0,1]. Values outside that\n\ |
5282 | 90 range are ignored.\n\ |
91 @end table\n\ | |
92 \n\ | |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
93 All other fields in @var{opts} are ignored. The outputs from @code{luinc}\n\ |
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
94 are the same as for @code{lu}.\n\ |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7243
diff
changeset
|
95 \n\ |
17281
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
15195
diff
changeset
|
96 Given the string argument @qcode{\"vector\"}, @code{luinc} returns the\n\ |
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
15195
diff
changeset
|
97 values of @var{p} @var{q} as vector values.\n\ |
19089
38937efbee21
Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents:
18851
diff
changeset
|
98 @seealso{sparse, lu, ilu, ichol}\n\ |
5642 | 99 @end deftypefn") |
5282 | 100 { |
101 int nargin = args.length (); | |
102 octave_value_list retval; | |
103 | |
104 if (nargin == 0) | |
5823 | 105 print_usage (); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7243
diff
changeset
|
106 else if (nargin < 2 || nargin > 3) |
5282 | 107 error ("luinc: incorrect number of arguments"); |
108 else | |
109 { | |
110 bool zero_level = false; | |
111 bool milu = false; | |
112 bool udiag = false; | |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7243
diff
changeset
|
113 Matrix thresh; |
5282 | 114 double droptol = -1.; |
13036
8afb81b32748
Initialise vecout variable and return permutation matrices instead of sparse matrices (bug #34185)
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents:
11593
diff
changeset
|
115 bool vecout = false; |
5282 | 116 |
117 if (args(1).is_string ()) | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
118 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
119 if (args(1).string_value () == "0") |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
120 zero_level = true; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
121 else |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
122 error ("luinc: unrecognized string argument"); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
123 } |
5282 | 124 else if (args(1).is_map ()) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
125 { |
11049
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
126 octave_scalar_map map = args(1).scalar_map_value (); |
5282 | 127 |
11049
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
128 if (! error_state) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
129 { |
11049
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
130 octave_value tmp; |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
131 |
11053
c33b7054f1f9
in recent Octave_map -> octave_scalar_map changes, use GETFIELD to access map elements, not CONTENTS
John W. Eaton <jwe@octave.org>
parents:
11049
diff
changeset
|
132 tmp = map.getfield ("droptol"); |
11049
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
133 if (tmp.is_defined ()) |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
134 droptol = tmp.double_value (); |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
135 |
11053
c33b7054f1f9
in recent Octave_map -> octave_scalar_map changes, use GETFIELD to access map elements, not CONTENTS
John W. Eaton <jwe@octave.org>
parents:
11049
diff
changeset
|
136 tmp = map.getfield ("milu"); |
11049
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
137 if (tmp.is_defined ()) |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
138 { |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
139 double val = tmp.double_value (); |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
140 |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
141 milu = (val == 0. ? false : true); |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
142 } |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
143 |
11053
c33b7054f1f9
in recent Octave_map -> octave_scalar_map changes, use GETFIELD to access map elements, not CONTENTS
John W. Eaton <jwe@octave.org>
parents:
11049
diff
changeset
|
144 tmp = map.getfield ("udiag"); |
11049
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
145 if (tmp.is_defined ()) |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
146 { |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
147 double val = tmp.double_value (); |
5282 | 148 |
11049
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
149 udiag = (val == 0. ? false : true); |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
150 } |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
151 |
11053
c33b7054f1f9
in recent Octave_map -> octave_scalar_map changes, use GETFIELD to access map elements, not CONTENTS
John W. Eaton <jwe@octave.org>
parents:
11049
diff
changeset
|
152 tmp = map.getfield ("thresh"); |
11049
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
153 if (tmp.is_defined ()) |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
154 { |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
155 thresh = tmp.matrix_value (); |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
156 |
20263
00cf2847355d
Deprecate Array::nelem() and Range::nelem() in favour of ::numel().
Carnë Draug <carandraug@octave.org>
parents:
20207
diff
changeset
|
157 if (thresh.numel () == 1) |
11049
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
158 { |
14854
5ae9f0f77635
maint: Use Octave coding conventions for coddling parenthis is DLD-FUNCTIONS directory
Rik <octave@nomad.inbox5.com>
parents:
14501
diff
changeset
|
159 thresh.resize (1,2); |
11049
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
160 thresh(1) = thresh(0); |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
161 } |
20263
00cf2847355d
Deprecate Array::nelem() and Range::nelem() in favour of ::numel().
Carnë Draug <carandraug@octave.org>
parents:
20207
diff
changeset
|
162 else if (thresh.numel () != 2) |
11049
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
163 { |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
164 error ("luinc: expecting 2-element vector for thresh"); |
11049
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
165 return retval; |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
166 } |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
167 } |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
168 } |
11049
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
169 else |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
170 { |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
171 error ("luinc: OPTS must be a scalar structure"); |
11049
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
172 return retval; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
173 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
174 } |
5282 | 175 else |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
176 droptol = args(1).double_value (); |
5282 | 177 |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7243
diff
changeset
|
178 if (nargin == 3) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
179 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
180 std::string tmp = args(2).string_value (); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7243
diff
changeset
|
181 |
20588
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
182 if (tmp.compare ("vector") == 0) |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
183 vecout = true; |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
184 else |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
185 error ("luinc: unrecognized string argument"); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
186 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7243
diff
changeset
|
187 |
17861
870f3e12e163
maint: Use phrase "FIXME:" for problem areas in code.
Rik <rik@octave.org>
parents:
17787
diff
changeset
|
188 // FIXME: Add code for zero-level factorization |
5282 | 189 if (zero_level) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
190 error ("luinc: zero-level factorization not implemented"); |
5282 | 191 |
20588
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
192 if (args(0).type_name () == "sparse matrix") |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
193 { |
20588
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
194 SparseMatrix sm = args(0).sparse_matrix_value (); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
195 octave_idx_type sm_nr = sm.rows (); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
196 octave_idx_type sm_nc = sm.cols (); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
197 ColumnVector Qinit (sm_nc); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
198 |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
199 for (octave_idx_type i = 0; i < sm_nc; i++) |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
200 Qinit (i) = i; |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
201 |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
202 switch (nargout) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
203 { |
20588
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
204 case 0: |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
205 case 1: |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
206 case 2: |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
207 { |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
208 SparseLU fact (sm, Qinit, thresh, false, true, droptol, |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
209 milu, udiag); |
5282 | 210 |
20588
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
211 SparseMatrix P = fact.Pr (); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
212 SparseMatrix L = P.transpose () * fact.L (); |
5282 | 213 |
20588
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
214 retval(1) |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
215 = octave_value (fact.U (), MatrixType (MatrixType::Upper)); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
216 |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
217 retval(0) |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
218 = octave_value (L, MatrixType (MatrixType::Permuted_Lower, |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
219 sm_nr, fact.row_perm ())); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
220 } |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
221 break; |
5282 | 222 |
20588
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
223 case 3: |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
224 { |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
225 SparseLU fact (sm, Qinit, thresh, false, true, droptol, |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
226 milu, udiag); |
5282 | 227 |
20588
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
228 if (vecout) |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
229 retval(2) = fact.Pr_vec (); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
230 else |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
231 retval(2) = fact.Pr_mat (); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
232 |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
233 retval(1) |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
234 = octave_value (fact.U (), MatrixType (MatrixType::Upper)); |
5282 | 235 |
20588
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
236 retval(0) |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
237 = octave_value (fact.L (), MatrixType (MatrixType::Lower)); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
238 } |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
239 break; |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
240 |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
241 case 4: |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
242 default: |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
243 { |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
244 SparseLU fact (sm, Qinit, thresh, false, false, droptol, |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
245 milu, udiag); |
5282 | 246 |
20588
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
247 if (vecout) |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
248 { |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
249 retval(3) = fact.Pc_vec (); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
250 retval(2) = fact.Pr_vec (); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
251 } |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
252 else |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
253 { |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
254 retval(3) = fact.Pc_mat (); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
255 retval(2) = fact.Pr_mat (); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
256 } |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
257 |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
258 retval(1) |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
259 = octave_value (fact.U (), MatrixType (MatrixType::Upper)); |
5282 | 260 |
20588
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
261 retval(0) |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
262 = octave_value (fact.L (), MatrixType (MatrixType::Lower)); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
263 } |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
264 break; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
265 } |
20588
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
266 } |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
267 else if (args(0).type_name () == "sparse complex matrix") |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
268 { |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
269 SparseComplexMatrix sm = |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
270 args(0).sparse_complex_matrix_value (); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
271 octave_idx_type sm_nr = sm.rows (); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
272 octave_idx_type sm_nc = sm.cols (); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
273 ColumnVector Qinit (sm_nc); |
5282 | 274 |
20588
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
275 for (octave_idx_type i = 0; i < sm_nc; i++) |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
276 Qinit (i) = i; |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
277 |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
278 switch (nargout) |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
279 { |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
280 case 0: |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
281 case 1: |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
282 case 2: |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
283 { |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
284 SparseComplexLU fact (sm, Qinit, thresh, false, true, |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
285 droptol, milu, udiag); |
5282 | 286 |
6027 | 287 |
20588
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
288 SparseMatrix P = fact.Pr (); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
289 SparseComplexMatrix L = P.transpose () * fact.L (); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
290 |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
291 retval(1) |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
292 = octave_value (fact.U (), MatrixType (MatrixType::Upper)); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
293 |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
294 retval(0) |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
295 = octave_value (L, MatrixType (MatrixType::Permuted_Lower, |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
296 sm_nr, fact.row_perm ())); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
297 } |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
298 break; |
5282 | 299 |
20588
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
300 case 3: |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
301 { |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
302 SparseComplexLU fact (sm, Qinit, thresh, false, true, |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
303 droptol, milu, udiag); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
304 |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
305 if (vecout) |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
306 retval(2) = fact.Pr_vec (); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
307 else |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
308 retval(2) = fact.Pr_mat (); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
309 |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
310 retval(1) |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
311 = octave_value (fact.U (), MatrixType (MatrixType::Upper)); |
5282 | 312 |
20588
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
313 retval(0) |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
314 = octave_value (fact.L (), MatrixType (MatrixType::Lower)); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
315 } |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
316 break; |
5282 | 317 |
20588
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
318 case 4: |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
319 default: |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
320 { |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
321 SparseComplexLU fact (sm, Qinit, thresh, false, false, |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
322 droptol, milu, udiag); |
5282 | 323 |
20588
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
324 if (vecout) |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
325 { |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
326 retval(3) = fact.Pc_vec (); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
327 retval(2) = fact.Pr_vec (); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
328 } |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
329 else |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
330 { |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
331 retval(3) = fact.Pc_mat (); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
332 retval(2) = fact.Pr_mat (); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
333 } |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
334 |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
335 retval(1) |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
336 = octave_value (fact.U (), MatrixType (MatrixType::Upper)); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
337 |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
338 retval(0) |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
339 = octave_value (fact.L (), MatrixType (MatrixType::Lower)); |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
340 } |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
341 break; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
342 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
343 } |
20588
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
344 else |
4bed806ee3d4
eliminate more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20263
diff
changeset
|
345 error ("luinc: matrix A must be sparse"); |
5282 | 346 } |
347 | |
348 return retval; | |
349 } | |
350 | |
351 /* | |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
352 %!testif HAVE_UMFPACK |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
353 %! a = sparse ([1,2,0,0;0,1,2,0;1e-14,0,3,0;0,0,0,1]); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
354 %! [l,u] = luinc (a, 1e-10); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
355 %! assert (l*u, sparse ([1,2,0,0;0,1,2,0;0,0,3,0;0,0,0,1]), 1e-10); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
356 %! opts.droptol = 1e-10; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
357 %! [l,u] = luinc (a, opts); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
358 %! assert (l*u, sparse ([1,2,0,0;0,1,2,0;0,0,3,0;0,0,0,1]), 1e-10); |
5610 | 359 |
7243 | 360 %!testif HAVE_UMFPACK |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
361 %! a = sparse ([1i,2,0,0;0,1,2,0;1e-14,0,3,0;0,0,0,1]); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
362 %! [l,u] = luinc (a, 1e-10); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
363 %! assert (l*u, sparse ([1i,2,0,0;0,1,2,0;0,0,3,0;0,0,0,1]), 1e-10); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
364 %! opts.droptol = 1e-10; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
365 %! [l,u] = luinc (a, opts); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
366 %! assert (l*u, sparse ([1i,2,0,0;0,1,2,0;0,0,3,0;0,0,0,1]), 1e-10); |
5610 | 367 */ |