5282
|
1 /* |
|
2 |
7017
|
3 Copyright (C) 2005, 2006, 2007 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 |
|
27 #include "defun-dld.h" |
|
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 |
|
40 DEFUN_DLD (luinc, args, nargout, |
|
41 "-*- texinfo -*-\n\ |
|
42 @deftypefn {Loadable Function} {[@var{l}, @var{u}, @var{p}, @var{q}] =} luinc (@var{a}, '0')\n\ |
|
43 @deftypefnx {Loadable Function} {[@var{l}, @var{u}, @var{p}, @var{q}] =} luinc (@var{a}, @var{droptol})\n\ |
|
44 @deftypefnx {Loadable Function} {[@var{l}, @var{u}, @var{p}, @var{q}] =} luinc (@var{a}, @var{opts})\n\ |
|
45 @cindex LU decomposition\n\ |
|
46 Produce the incomplete LU factorization of the sparse matrix @var{a}.\n\ |
|
47 Two types of incomplete factorization are possible, and the type\n\ |
|
48 is determined by the second argument to @dfn{luinc}.\n\ |
|
49 \n\ |
|
50 Called with a second argument of '0', the zero-level incomplete\n\ |
|
51 LU factorization is produced. This creates a factorization of @var{a}\n\ |
|
52 where the position of the non-zero arguments correspond to the same\n\ |
|
53 positions as in the matrix @var{a}.\n\ |
|
54 \n\ |
|
55 Alternatively, the fill-in of the incomplete LU factorization can\n\ |
|
56 be controlled through the variable @var{droptol} or the structure\n\ |
|
57 @var{opts}. The UMFPACK multifrontal factorization code by Tim A.\n\ |
7001
|
58 Davis is used for the incomplete LU factorization, (availability\n\ |
5282
|
59 @url{http://www.cise.ufl.edu/research/sparse/umfpack/})\n\ |
|
60 \n\ |
|
61 @var{droptol} determines the values below which the values in the LU\n\ |
|
62 factorization are dropped and replaced by zero. It must be a positive\n\ |
|
63 scalar, and any values in the factorization whose absolute value are\n\ |
|
64 less than this value are dropped, expect if leaving them increase the\n\ |
|
65 sparsity of the matrix. Setting @var{droptol} to zero results in a\n\ |
|
66 complete LU factorization which is the default.\n\ |
|
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\ |
|
72 The drop tolerance as above. If @var{opts} only contains @code{droptol}\n\ |
|
73 then this is equivalent to using the variable @var{droptol}.\n\ |
|
74 \n\ |
|
75 @item milu\n\ |
|
76 A logical variable flagging whether to use the modified incomplete LU\n\ |
|
77 factorization. In the case that @code{milu} is true, the dropped values\n\ |
7001
|
78 are subtracted from the diagonal of the matrix U of the factorization.\n\ |
5282
|
79 The default is @code{false}.\n\ |
|
80 \n\ |
|
81 @item udiag\n\ |
|
82 A logical variable that flags whether zero elements on the diagonal of U\n\ |
|
83 should be replaced with @var{droptol} to attempt to avoid singular\n\ |
|
84 factors. The default is @code{false}.\n\ |
|
85 \n\ |
|
86 @item thresh\n\ |
|
87 Defines the pivot threshold in the interval [0,1]. Values outside that\n\ |
|
88 range are ignored.\n\ |
|
89 @end table\n\ |
|
90 \n\ |
|
91 All other fields in @var{opts} are ignored. The outputs from @dfn{luinc}\n\ |
|
92 are the same as for @dfn{lu}.\n\ |
5642
|
93 @seealso{sparse, lu, cholinc}\n\ |
|
94 @end deftypefn") |
5282
|
95 { |
|
96 int nargin = args.length (); |
|
97 octave_value_list retval; |
|
98 |
|
99 if (nargin == 0) |
5823
|
100 print_usage (); |
5282
|
101 else if (nargin != 2) |
|
102 error ("luinc: incorrect number of arguments"); |
|
103 else |
|
104 { |
|
105 bool zero_level = false; |
|
106 bool milu = false; |
|
107 bool udiag = false; |
|
108 bool thresh = -1; |
|
109 double droptol = -1.; |
|
110 |
|
111 if (args(1).is_string ()) |
|
112 { |
|
113 if (args(1).string_value () == "0") |
|
114 zero_level = true; |
|
115 else |
|
116 error ("luinc: unrecognized string argument"); |
|
117 } |
|
118 else if (args(1).is_map ()) |
|
119 { |
|
120 Octave_map map = args(1).map_value (); |
|
121 |
|
122 if (map.contains ("droptol")) |
|
123 droptol = map.contents ("droptol")(0).double_value (); |
|
124 |
|
125 if (map.contains ("milu")) |
|
126 { |
|
127 double tmp = map.contents ("milu")(0).double_value (); |
|
128 |
|
129 milu = (tmp == 0. ? false : true); |
|
130 } |
|
131 |
|
132 if (map.contains ("udiag")) |
|
133 { |
|
134 double tmp = map.contents ("udiag")(0).double_value (); |
|
135 |
|
136 udiag = (tmp == 0. ? false : true); |
|
137 } |
|
138 |
|
139 if (map.contains ("thresh")) |
|
140 thresh = map.contents ("thresh")(0).double_value (); |
|
141 } |
|
142 else |
|
143 droptol = args(1).double_value (); |
|
144 |
5775
|
145 // FIXME Add code for zero-level factorization |
5282
|
146 if (zero_level) |
|
147 error ("luinc: zero-level factorization not implemented"); |
|
148 |
|
149 if (!error_state) |
|
150 { |
5322
|
151 if (args(0).type_name () == "sparse matrix") |
5282
|
152 { |
|
153 SparseMatrix sm = args(0).sparse_matrix_value (); |
5322
|
154 octave_idx_type sm_nr = sm.rows (); |
5282
|
155 octave_idx_type sm_nc = sm.cols (); |
|
156 ColumnVector Qinit (sm_nc); |
|
157 |
|
158 for (octave_idx_type i = 0; i < sm_nc; i++) |
|
159 Qinit (i) = i; |
|
160 |
|
161 if (! error_state) |
|
162 { |
|
163 switch (nargout) |
|
164 { |
|
165 case 0: |
|
166 case 1: |
|
167 case 2: |
|
168 { |
|
169 SparseLU fact (sm, Qinit, thresh, true, droptol, |
|
170 milu, udiag); |
|
171 |
6027
|
172 if (! error_state) |
|
173 { |
|
174 SparseMatrix P = fact.Pr (); |
|
175 SparseMatrix L = P.transpose () * fact.L (); |
|
176 retval(1) = octave_value (fact.U (), |
|
177 MatrixType (MatrixType::Upper)); |
|
178 retval(0) = octave_value (L, MatrixType |
|
179 (MatrixType::Permuted_Lower, |
|
180 sm_nr, fact.row_perm ())); |
|
181 } |
5282
|
182 } |
|
183 break; |
|
184 |
|
185 case 3: |
|
186 { |
|
187 SparseLU fact (sm, Qinit, thresh, true, droptol, |
|
188 milu, udiag); |
|
189 |
6027
|
190 if (! error_state) |
|
191 { |
|
192 retval(2) = fact.Pr (); |
|
193 retval(1) = octave_value (fact.U (), |
|
194 MatrixType (MatrixType::Upper)); |
|
195 retval(0) = octave_value (fact.L (), |
|
196 MatrixType (MatrixType::Lower)); |
|
197 } |
5282
|
198 } |
|
199 break; |
|
200 |
|
201 case 4: |
|
202 default: |
|
203 { |
|
204 SparseLU fact (sm, Qinit, thresh, false, droptol, |
|
205 milu, udiag); |
|
206 |
6027
|
207 if (! error_state) |
|
208 { |
|
209 retval(3) = fact.Pc (); |
|
210 retval(2) = fact.Pr (); |
|
211 retval(1) = octave_value (fact.U (), |
|
212 MatrixType (MatrixType::Upper)); |
|
213 retval(0) = octave_value (fact.L (), |
|
214 MatrixType (MatrixType::Lower)); |
|
215 } |
5282
|
216 } |
|
217 break; |
|
218 } |
|
219 } |
|
220 } |
5322
|
221 else if (args(0).type_name () == "sparse complex matrix") |
5282
|
222 { |
|
223 SparseComplexMatrix sm = |
|
224 args(0).sparse_complex_matrix_value (); |
5322
|
225 octave_idx_type sm_nr = sm.rows (); |
5282
|
226 octave_idx_type sm_nc = sm.cols (); |
|
227 ColumnVector Qinit (sm_nc); |
|
228 |
|
229 for (octave_idx_type i = 0; i < sm_nc; i++) |
|
230 Qinit (i) = i; |
|
231 |
|
232 if (! error_state) |
|
233 { |
|
234 switch (nargout) |
|
235 { |
|
236 case 0: |
|
237 case 1: |
|
238 case 2: |
|
239 { |
|
240 SparseComplexLU fact (sm, Qinit, thresh, true, |
|
241 droptol, milu, udiag); |
|
242 |
6027
|
243 |
|
244 if (! error_state) |
|
245 { |
|
246 SparseMatrix P = fact.Pr (); |
|
247 SparseComplexMatrix L = P.transpose () * fact.L (); |
|
248 retval(1) = octave_value (fact.U (), |
|
249 MatrixType (MatrixType::Upper)); |
|
250 retval(0) = octave_value (L, MatrixType |
|
251 (MatrixType::Permuted_Lower, |
|
252 sm_nr, fact.row_perm ())); |
|
253 } |
5282
|
254 } |
|
255 break; |
|
256 |
|
257 case 3: |
|
258 { |
|
259 SparseComplexLU fact (sm, Qinit, thresh, true, |
|
260 droptol, milu, udiag); |
|
261 |
6027
|
262 if (! error_state) |
|
263 { |
|
264 retval(2) = fact.Pr (); |
|
265 retval(1) = octave_value (fact.U (), |
|
266 MatrixType (MatrixType::Upper)); |
|
267 retval(0) = octave_value (fact.L (), |
|
268 MatrixType (MatrixType::Lower)); |
|
269 } |
5282
|
270 } |
|
271 break; |
|
272 |
|
273 case 4: |
|
274 default: |
|
275 { |
|
276 SparseComplexLU fact (sm, Qinit, thresh, false, |
|
277 droptol, milu, udiag); |
|
278 |
6027
|
279 if (! error_state) |
|
280 { |
|
281 retval(3) = fact.Pc (); |
|
282 retval(2) = fact.Pr (); |
|
283 retval(1) = octave_value (fact.U (), |
|
284 MatrixType (MatrixType::Upper)); |
|
285 retval(0) = octave_value (fact.L (), |
|
286 MatrixType (MatrixType::Lower)); |
|
287 } |
5282
|
288 } |
|
289 break; |
|
290 } |
|
291 } |
|
292 } |
|
293 else |
|
294 error ("luinc: first argument must be sparse"); |
|
295 } |
|
296 } |
|
297 |
|
298 return retval; |
|
299 } |
|
300 |
|
301 /* |
5610
|
302 |
|
303 %!test |
|
304 %! a=sparse([1,2,0,0;0,1,2,0;1e-14,0,3,0;0,0,0,1]); |
|
305 %! [l,u]=luinc(a,1e-10); |
|
306 %! assert(l*u, sparse([1,2,0,0;0,1,2,0;0,0,3,0;0,0,0,1]),1e-10); |
|
307 %! opts.droptol=1e-10; |
|
308 %! [l,u]=luinc(a,opts); |
|
309 %! assert(l*u, sparse([1,2,0,0;0,1,2,0;0,0,3,0;0,0,0,1]),1e-10); |
|
310 |
|
311 %!test |
|
312 %! a=sparse([1i,2,0,0;0,1,2,0;1e-14,0,3,0;0,0,0,1]); |
|
313 %! [l,u]=luinc(a,1e-10); |
|
314 %! assert(l*u, sparse([1i,2,0,0;0,1,2,0;0,0,3,0;0,0,0,1]),1e-10); |
|
315 %! opts.droptol=1e-10; |
|
316 %! [l,u]=luinc(a,opts); |
|
317 %! assert(l*u, sparse([1i,2,0,0;0,1,2,0;0,0,3,0;0,0,0,1]),1e-10); |
|
318 |
|
319 */ |
|
320 |
|
321 /* |
5282
|
322 ;;; Local Variables: *** |
|
323 ;;; mode: C++ *** |
|
324 ;;; End: *** |
|
325 */ |