Mercurial > hg > octave-nkf
annotate libinterp/corefcn/luinc.cc @ 19840:c5270263d466 gui-release
close gui-release branch
author | John W. Eaton <jwe@octave.org> |
---|---|
date | Fri, 30 Jan 2015 17:41:50 -0500 |
parents | 870f3e12e163 |
children | 6113e0c6920b |
rev | line source |
---|---|
5282 | 1 /* |
2 | |
17744
d63878346099
maint: Update copyright notices for release.
John W. Eaton <jwe@octave.org>
parents:
17281
diff
changeset
|
3 Copyright (C) 2005-2013 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 | |
15039
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14854
diff
changeset
|
40 DEFUN (luinc, args, nargout, |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
41 "-*- texinfo -*-\n\ |
15039
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14854
diff
changeset
|
42 @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
|
43 @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
|
44 @deftypefnx {Built-in Function} {[@var{L}, @var{U}, @var{P}, @var{Q}] =} luinc (@var{A}, @var{opts})\n\ |
5282 | 45 @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
|
46 Produce the incomplete LU@tie{}factorization of the sparse matrix @var{A}.\n\ |
5282 | 47 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
|
48 is determined by the second argument to @code{luinc}.\n\ |
5282 | 49 \n\ |
17281
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
15195
diff
changeset
|
50 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
|
51 LU@tie{}factorization is produced. This creates a factorization of @var{A}\n\ |
5282 | 52 where the position of the non-zero 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
|
53 positions as in the matrix @var{A}.\n\ |
5282 | 54 \n\ |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
55 Alternatively, the fill-in of the incomplete LU@tie{}factorization can\n\ |
5282 | 56 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
|
57 @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
|
58 Davis is used for the incomplete LU@tie{}factorization, (availability\n\ |
5282 | 59 @url{http://www.cise.ufl.edu/research/sparse/umfpack/})\n\ |
60 \n\ | |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
61 @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
|
62 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
|
63 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
|
64 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
|
65 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
|
66 LU@tie{}factorization which is the default.\n\ |
5282 | 67 \n\ |
68 @var{opts} is a structure containing one or more of the fields\n\ | |
69 \n\ | |
70 @table @code\n\ | |
71 @item droptol\n\ | |
9064
7c02ec148a3c
Check grammar on all .cc files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
72 The drop tolerance as above. If @var{opts} only contains @code{droptol}\n\ |
5282 | 73 then this is equivalent to using the variable @var{droptol}.\n\ |
74 \n\ | |
75 @item milu\n\ | |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
76 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
|
77 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
|
78 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
|
79 factorization. The default is @code{false}.\n\ |
5282 | 80 \n\ |
81 @item udiag\n\ | |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
82 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
|
83 @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
|
84 factors. The default is @code{false}.\n\ |
5282 | 85 \n\ |
86 @item thresh\n\ | |
9064
7c02ec148a3c
Check grammar on all .cc files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
87 Defines the pivot threshold in the interval [0,1]. Values outside that\n\ |
5282 | 88 range are ignored.\n\ |
89 @end table\n\ | |
90 \n\ | |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
91 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
|
92 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
|
93 \n\ |
17281
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
15195
diff
changeset
|
94 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
|
95 values of @var{p} @var{q} as vector values.\n\ |
8286
6f2d95255911
fix @seealso references to point to existing anchors
Thorsten Meyer <thorsten.meyier@gmx.de>
parents:
7515
diff
changeset
|
96 @seealso{sparse, lu}\n\ |
5642 | 97 @end deftypefn") |
5282 | 98 { |
99 int nargin = args.length (); | |
100 octave_value_list retval; | |
101 | |
102 if (nargin == 0) | |
5823 | 103 print_usage (); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7243
diff
changeset
|
104 else if (nargin < 2 || nargin > 3) |
5282 | 105 error ("luinc: incorrect number of arguments"); |
106 else | |
107 { | |
108 bool zero_level = false; | |
109 bool milu = false; | |
110 bool udiag = false; | |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7243
diff
changeset
|
111 Matrix thresh; |
5282 | 112 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
|
113 bool vecout = false; |
5282 | 114 |
115 if (args(1).is_string ()) | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
116 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
117 if (args(1).string_value () == "0") |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
118 zero_level = true; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
119 else |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
120 error ("luinc: unrecognized string argument"); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
121 } |
5282 | 122 else if (args(1).is_map ()) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
123 { |
11049
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
124 octave_scalar_map map = args(1).scalar_map_value (); |
5282 | 125 |
11049
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
126 if (! error_state) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
127 { |
11049
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
128 octave_value tmp; |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
129 |
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
|
130 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
|
131 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
|
132 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
|
133 |
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
|
134 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
|
135 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
|
136 { |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
137 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
|
138 |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
139 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
|
140 } |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
141 |
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
|
142 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
|
143 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
|
144 { |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
145 double val = tmp.double_value (); |
5282 | 146 |
11049
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
147 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
|
148 } |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
149 |
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
|
150 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
|
151 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
|
152 { |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
153 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
|
154 |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
155 if (thresh.nelem () == 1) |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
156 { |
14854
5ae9f0f77635
maint: Use Octave coding conventions for coddling parenthis is DLD-FUNCTIONS directory
Rik <octave@nomad.inbox5.com>
parents:
14501
diff
changeset
|
157 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
|
158 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
|
159 } |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
160 else if (thresh.nelem () != 2) |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
161 { |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
162 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
|
163 return retval; |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
164 } |
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
165 } |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
166 } |
11049
b0a9450d81c6
luinc.cc: use octave_scalar_map instead of Octave_map
John W. Eaton <jwe@octave.org>
parents:
10840
diff
changeset
|
167 else |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
168 { |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
169 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
|
170 return retval; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
171 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
172 } |
5282 | 173 else |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
174 droptol = args(1).double_value (); |
5282 | 175 |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7243
diff
changeset
|
176 if (nargin == 3) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
177 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
178 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
|
179 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
180 if (! error_state ) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
181 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
182 if (tmp.compare ("vector") == 0) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
183 vecout = true; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
184 else |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
185 error ("luinc: unrecognized string argument"); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
186 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
187 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7243
diff
changeset
|
188 |
17861
870f3e12e163
maint: Use phrase "FIXME:" for problem areas in code.
Rik <rik@octave.org>
parents:
17787
diff
changeset
|
189 // FIXME: Add code for zero-level factorization |
5282 | 190 if (zero_level) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
191 error ("luinc: zero-level factorization not implemented"); |
5282 | 192 |
193 if (!error_state) | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
194 { |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11553
diff
changeset
|
195 if (args(0).type_name () == "sparse matrix") |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
196 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
197 SparseMatrix sm = args(0).sparse_matrix_value (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
198 octave_idx_type sm_nr = sm.rows (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
199 octave_idx_type sm_nc = sm.cols (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
200 ColumnVector Qinit (sm_nc); |
5282 | 201 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
202 for (octave_idx_type i = 0; i < sm_nc; i++) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
203 Qinit (i) = i; |
5282 | 204 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
205 if (! error_state) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
206 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
207 switch (nargout) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
208 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
209 case 0: |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
210 case 1: |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
211 case 2: |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
212 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
213 SparseLU fact (sm, Qinit, thresh, false, true, droptol, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
214 milu, udiag); |
5282 | 215 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
216 if (! error_state) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
217 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
218 SparseMatrix P = fact.Pr (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
219 SparseMatrix L = P.transpose () * fact.L (); |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
220 retval(1) |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
221 = octave_value (fact.U (), |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
222 MatrixType (MatrixType::Upper)); |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
223 retval(0) |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
224 = octave_value (L, MatrixType |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
225 (MatrixType::Permuted_Lower, |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
226 sm_nr, fact.row_perm ())); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
227 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
228 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
229 break; |
5282 | 230 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
231 case 3: |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
232 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
233 SparseLU fact (sm, Qinit, thresh, false, true, droptol, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
234 milu, udiag); |
5282 | 235 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
236 if (! error_state) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
237 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
238 if (vecout) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
239 retval(2) = fact.Pr_vec (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
240 else |
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
|
241 retval(2) = fact.Pr_mat (); |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
242 retval(1) |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
243 = octave_value (fact.U (), |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
244 MatrixType (MatrixType::Upper)); |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
245 retval(0) |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
246 = octave_value (fact.L (), |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
247 MatrixType (MatrixType::Lower)); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
248 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
249 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
250 break; |
5282 | 251 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
252 case 4: |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
253 default: |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
254 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
255 SparseLU fact (sm, Qinit, thresh, false, false, droptol, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
256 milu, udiag); |
5282 | 257 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
258 if (! error_state) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
259 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
260 if (vecout) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
261 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
262 retval(3) = fact.Pc_vec (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
263 retval(2) = fact.Pr_vec (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
264 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
265 else |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
266 { |
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
|
267 retval(3) = fact.Pc_mat (); |
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
|
268 retval(2) = fact.Pr_mat (); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
269 } |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
270 retval(1) |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
271 = octave_value (fact.U (), |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
272 MatrixType (MatrixType::Upper)); |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
273 retval(0) |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
274 = octave_value (fact.L (), |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
275 MatrixType (MatrixType::Lower)); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
276 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
277 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
278 break; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
279 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
280 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
281 } |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11553
diff
changeset
|
282 else if (args(0).type_name () == "sparse complex matrix") |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
283 { |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11553
diff
changeset
|
284 SparseComplexMatrix sm = |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
285 args(0).sparse_complex_matrix_value (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
286 octave_idx_type sm_nr = sm.rows (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
287 octave_idx_type sm_nc = sm.cols (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
288 ColumnVector Qinit (sm_nc); |
5282 | 289 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
290 for (octave_idx_type i = 0; i < sm_nc; i++) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
291 Qinit (i) = i; |
5282 | 292 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
293 if (! error_state) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
294 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
295 switch (nargout) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
296 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
297 case 0: |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
298 case 1: |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
299 case 2: |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
300 { |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11553
diff
changeset
|
301 SparseComplexLU fact (sm, Qinit, thresh, false, true, |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
302 droptol, milu, udiag); |
5282 | 303 |
6027 | 304 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
305 if (! error_state) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
306 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
307 SparseMatrix P = fact.Pr (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
308 SparseComplexMatrix L = P.transpose () * fact.L (); |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
309 retval(1) |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
310 = octave_value (fact.U (), |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
311 MatrixType (MatrixType::Upper)); |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
312 retval(0) |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
313 = octave_value (L, MatrixType |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
314 (MatrixType::Permuted_Lower, |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
315 sm_nr, fact.row_perm ())); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
316 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
317 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
318 break; |
5282 | 319 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
320 case 3: |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
321 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
322 SparseComplexLU fact (sm, Qinit, thresh, false, true, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
323 droptol, milu, udiag); |
5282 | 324 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
325 if (! error_state) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
326 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
327 if (vecout) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
328 retval(2) = fact.Pr_vec (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
329 else |
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
|
330 retval(2) = fact.Pr_mat (); |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
331 retval(1) |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
332 = octave_value (fact.U (), |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
333 MatrixType (MatrixType::Upper)); |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
334 retval(0) |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
335 = octave_value (fact.L (), |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
336 MatrixType (MatrixType::Lower)); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
337 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
338 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
339 break; |
5282 | 340 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
341 case 4: |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
342 default: |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
343 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
344 SparseComplexLU fact (sm, Qinit, thresh, false, false, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
345 droptol, milu, udiag); |
5282 | 346 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
347 if (! error_state) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
348 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
349 if (vecout) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
350 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
351 retval(3) = fact.Pc_vec (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
352 retval(2) = fact.Pr_vec (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
353 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
354 else |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
355 { |
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
|
356 retval(3) = fact.Pc_mat (); |
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
|
357 retval(2) = fact.Pr_mat (); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
358 } |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
359 retval(1) |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
360 = octave_value (fact.U (), |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
361 MatrixType (MatrixType::Upper)); |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
362 retval(0) |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
363 = octave_value (fact.L (), |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
364 MatrixType (MatrixType::Lower)); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
365 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
366 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
367 break; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
368 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
369 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
370 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
371 else |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
372 error ("luinc: matrix A must be sparse"); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9064
diff
changeset
|
373 } |
5282 | 374 } |
375 | |
376 return retval; | |
377 } | |
378 | |
379 /* | |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
380 %!testif HAVE_UMFPACK |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
381 %! 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
|
382 %! [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
|
383 %! 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
|
384 %! opts.droptol = 1e-10; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
385 %! [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
|
386 %! assert (l*u, sparse ([1,2,0,0;0,1,2,0;0,0,3,0;0,0,0,1]), 1e-10); |
5610 | 387 |
7243 | 388 %!testif HAVE_UMFPACK |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
389 %! 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
|
390 %! [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
|
391 %! 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
|
392 %! opts.droptol = 1e-10; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
393 %! [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
|
394 %! assert (l*u, sparse ([1i,2,0,0;0,1,2,0;0,0,3,0;0,0,0,1]), 1e-10); |
5610 | 395 */ |