5282
|
1 /* |
|
2 |
|
3 Copyright (C) 2005 David Bateman |
|
4 |
|
5 Octave is free software; you can redistribute it and/or modify it |
|
6 under the terms of the GNU General Public License as published by the |
|
7 Free Software Foundation; either version 2, or (at your option) any |
|
8 later version. |
|
9 |
|
10 Octave is distributed in the hope that it will be useful, but WITHOUT |
|
11 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or |
|
12 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License |
|
13 for more details. |
|
14 |
|
15 You should have received a copy of the GNU General Public License |
5307
|
16 along with this program; see the file COPYING. If not, write to the |
|
17 Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, |
|
18 Boston, MA 02110-1301, USA. |
5282
|
19 |
|
20 */ |
|
21 |
|
22 #ifdef HAVE_CONFIG_H |
|
23 #include <config.h> |
|
24 #endif |
|
25 |
|
26 #include "defun-dld.h" |
|
27 #include "error.h" |
|
28 #include "gripes.h" |
|
29 #include "oct-obj.h" |
|
30 #include "utils.h" |
|
31 #include "oct-map.h" |
|
32 |
5322
|
33 #include "SparseType.h" |
5282
|
34 #include "SparseCmplxLU.h" |
|
35 #include "SparsedbleLU.h" |
|
36 #include "ov-re-sparse.h" |
|
37 #include "ov-cx-sparse.h" |
|
38 |
|
39 DEFUN_DLD (luinc, args, nargout, |
|
40 "-*- texinfo -*-\n\ |
|
41 @deftypefn {Loadable Function} {[@var{l}, @var{u}, @var{p}, @var{q}] =} luinc (@var{a}, '0')\n\ |
|
42 @deftypefnx {Loadable Function} {[@var{l}, @var{u}, @var{p}, @var{q}] =} luinc (@var{a}, @var{droptol})\n\ |
|
43 @deftypefnx {Loadable Function} {[@var{l}, @var{u}, @var{p}, @var{q}] =} luinc (@var{a}, @var{opts})\n\ |
|
44 @cindex LU decomposition\n\ |
|
45 Produce the incomplete LU factorization of the sparse matrix @var{a}.\n\ |
|
46 Two types of incomplete factorization are possible, and the type\n\ |
|
47 is determined by the second argument to @dfn{luinc}.\n\ |
|
48 \n\ |
|
49 Called with a second argument of '0', the zero-level incomplete\n\ |
|
50 LU factorization is produced. This creates a factorization of @var{a}\n\ |
|
51 where the position of the non-zero arguments correspond to the same\n\ |
|
52 positions as in the matrix @var{a}.\n\ |
|
53 \n\ |
|
54 Alternatively, the fill-in of the incomplete LU factorization can\n\ |
|
55 be controlled through the variable @var{droptol} or the structure\n\ |
|
56 @var{opts}. The UMFPACK multifrontal factorization code by Tim A.\n\ |
|
57 Davis is used for the incomplete LU factorication, (availability\n\ |
|
58 @url{http://www.cise.ufl.edu/research/sparse/umfpack/})\n\ |
|
59 \n\ |
|
60 @var{droptol} determines the values below which the values in the LU\n\ |
|
61 factorization are dropped and replaced by zero. It must be a positive\n\ |
|
62 scalar, and any values in the factorization whose absolute value are\n\ |
|
63 less than this value are dropped, expect if leaving them increase the\n\ |
|
64 sparsity of the matrix. Setting @var{droptol} to zero results in a\n\ |
|
65 complete LU factorization which is the default.\n\ |
|
66 \n\ |
|
67 @var{opts} is a structure containing one or more of the fields\n\ |
|
68 \n\ |
|
69 @table @code\n\ |
|
70 @item droptol\n\ |
|
71 The drop tolerance as above. If @var{opts} only contains @code{droptol}\n\ |
|
72 then this is equivalent to using the variable @var{droptol}.\n\ |
|
73 \n\ |
|
74 @item milu\n\ |
|
75 A logical variable flagging whether to use the modified incomplete LU\n\ |
|
76 factorization. In the case that @code{milu} is true, the dropped values\n\ |
|
77 are subtract from the diagonal of the matrix U of the factorization.\n\ |
|
78 The default is @code{false}.\n\ |
|
79 \n\ |
|
80 @item udiag\n\ |
|
81 A logical variable that flags whether zero elements on the diagonal of U\n\ |
|
82 should be replaced with @var{droptol} to attempt to avoid singular\n\ |
|
83 factors. The default is @code{false}.\n\ |
|
84 \n\ |
|
85 @item thresh\n\ |
|
86 Defines the pivot threshold in the interval [0,1]. Values outside that\n\ |
|
87 range are ignored.\n\ |
|
88 @end table\n\ |
|
89 \n\ |
|
90 All other fields in @var{opts} are ignored. The outputs from @dfn{luinc}\n\ |
|
91 are the same as for @dfn{lu}.\n\ |
|
92 @end deftypefn\n\ |
|
93 @seealso{sparse, lu, cholinc}") |
|
94 { |
|
95 int nargin = args.length (); |
|
96 octave_value_list retval; |
|
97 |
|
98 if (nargin == 0) |
|
99 print_usage ("luinc"); |
|
100 else if (nargin != 2) |
|
101 error ("luinc: incorrect number of arguments"); |
|
102 else |
|
103 { |
|
104 bool zero_level = false; |
|
105 bool milu = false; |
|
106 bool udiag = false; |
|
107 bool thresh = -1; |
|
108 double droptol = -1.; |
|
109 |
|
110 if (args(1).is_string ()) |
|
111 { |
|
112 if (args(1).string_value () == "0") |
|
113 zero_level = true; |
|
114 else |
|
115 error ("luinc: unrecognized string argument"); |
|
116 } |
|
117 else if (args(1).is_map ()) |
|
118 { |
|
119 Octave_map map = args(1).map_value (); |
|
120 |
|
121 if (map.contains ("droptol")) |
|
122 droptol = map.contents ("droptol")(0).double_value (); |
|
123 |
|
124 if (map.contains ("milu")) |
|
125 { |
|
126 double tmp = map.contents ("milu")(0).double_value (); |
|
127 |
|
128 milu = (tmp == 0. ? false : true); |
|
129 } |
|
130 |
|
131 if (map.contains ("udiag")) |
|
132 { |
|
133 double tmp = map.contents ("udiag")(0).double_value (); |
|
134 |
|
135 udiag = (tmp == 0. ? false : true); |
|
136 } |
|
137 |
|
138 if (map.contains ("thresh")) |
|
139 thresh = map.contents ("thresh")(0).double_value (); |
|
140 } |
|
141 else |
|
142 droptol = args(1).double_value (); |
|
143 |
|
144 // XXX FIXME XXX Add code for zero-level factorization |
|
145 if (zero_level) |
|
146 error ("luinc: zero-level factorization not implemented"); |
|
147 |
|
148 if (!error_state) |
|
149 { |
5322
|
150 if (args(0).type_name () == "sparse matrix") |
5282
|
151 { |
|
152 SparseMatrix sm = args(0).sparse_matrix_value (); |
5322
|
153 octave_idx_type sm_nr = sm.rows (); |
5282
|
154 octave_idx_type sm_nc = sm.cols (); |
|
155 ColumnVector Qinit (sm_nc); |
|
156 |
|
157 for (octave_idx_type i = 0; i < sm_nc; i++) |
|
158 Qinit (i) = i; |
|
159 |
|
160 if (! error_state) |
|
161 { |
|
162 switch (nargout) |
|
163 { |
|
164 case 0: |
|
165 case 1: |
|
166 case 2: |
|
167 { |
|
168 SparseLU fact (sm, Qinit, thresh, true, droptol, |
|
169 milu, udiag); |
|
170 |
|
171 SparseMatrix P = fact.Pr (); |
|
172 SparseMatrix L = P.transpose () * fact.L (); |
5322
|
173 retval(1) = octave_value (fact.U (), |
|
174 SparseType (SparseType::Upper)); |
|
175 retval(0) = octave_value (L, SparseType |
|
176 (SparseType::Permuted_Lower, |
|
177 sm_nr, fact.row_perm ())); |
5282
|
178 } |
|
179 break; |
|
180 |
|
181 case 3: |
|
182 { |
|
183 SparseLU fact (sm, Qinit, thresh, true, droptol, |
|
184 milu, udiag); |
|
185 |
|
186 retval(2) = fact.Pr (); |
5322
|
187 retval(1) = octave_value (fact.U (), |
|
188 SparseType (SparseType::Upper)); |
|
189 retval(0) = octave_value (fact.L (), |
|
190 SparseType (SparseType::Lower)); |
5282
|
191 } |
|
192 break; |
|
193 |
|
194 case 4: |
|
195 default: |
|
196 { |
|
197 SparseLU fact (sm, Qinit, thresh, false, droptol, |
|
198 milu, udiag); |
|
199 |
|
200 retval(3) = fact.Pc (); |
|
201 retval(2) = fact.Pr (); |
5322
|
202 retval(1) = octave_value (fact.U (), |
|
203 SparseType (SparseType::Upper)); |
|
204 retval(0) = octave_value (fact.L (), |
|
205 SparseType (SparseType::Lower)); |
5282
|
206 } |
|
207 break; |
|
208 } |
|
209 } |
|
210 } |
5322
|
211 else if (args(0).type_name () == "sparse complex matrix") |
5282
|
212 { |
|
213 SparseComplexMatrix sm = |
|
214 args(0).sparse_complex_matrix_value (); |
5322
|
215 octave_idx_type sm_nr = sm.rows (); |
5282
|
216 octave_idx_type sm_nc = sm.cols (); |
|
217 ColumnVector Qinit (sm_nc); |
|
218 |
|
219 for (octave_idx_type i = 0; i < sm_nc; i++) |
|
220 Qinit (i) = i; |
|
221 |
|
222 if (! error_state) |
|
223 { |
|
224 switch (nargout) |
|
225 { |
|
226 case 0: |
|
227 case 1: |
|
228 case 2: |
|
229 { |
|
230 SparseComplexLU fact (sm, Qinit, thresh, true, |
|
231 droptol, milu, udiag); |
|
232 |
|
233 SparseMatrix P = fact.Pr (); |
|
234 SparseComplexMatrix L = P.transpose () * fact.L (); |
5322
|
235 retval(1) = octave_value (fact.U (), |
|
236 SparseType (SparseType::Upper)); |
|
237 retval(0) = octave_value (L, SparseType |
|
238 (SparseType::Permuted_Lower, |
|
239 sm_nr, fact.row_perm ())); |
5282
|
240 } |
|
241 break; |
|
242 |
|
243 case 3: |
|
244 { |
|
245 SparseComplexLU fact (sm, Qinit, thresh, true, |
|
246 droptol, milu, udiag); |
|
247 |
|
248 retval(2) = fact.Pr (); |
5322
|
249 retval(1) = octave_value (fact.U (), |
|
250 SparseType (SparseType::Upper)); |
|
251 retval(0) = octave_value (fact.L (), |
|
252 SparseType (SparseType::Lower)); |
5282
|
253 } |
|
254 break; |
|
255 |
|
256 case 4: |
|
257 default: |
|
258 { |
|
259 SparseComplexLU fact (sm, Qinit, thresh, false, |
|
260 droptol, milu, udiag); |
|
261 |
|
262 retval(3) = fact.Pc (); |
|
263 retval(2) = fact.Pr (); |
5322
|
264 retval(1) = octave_value (fact.U (), |
|
265 SparseType (SparseType::Upper)); |
|
266 retval(0) = octave_value (fact.L (), |
|
267 SparseType (SparseType::Lower)); |
5282
|
268 } |
|
269 break; |
|
270 } |
|
271 } |
|
272 } |
|
273 else |
|
274 error ("luinc: first argument must be sparse"); |
|
275 } |
|
276 } |
|
277 |
|
278 return retval; |
|
279 } |
|
280 |
|
281 /* |
5610
|
282 |
|
283 %!test |
|
284 %! a=sparse([1,2,0,0;0,1,2,0;1e-14,0,3,0;0,0,0,1]); |
|
285 %! [l,u]=luinc(a,1e-10); |
|
286 %! assert(l*u, sparse([1,2,0,0;0,1,2,0;0,0,3,0;0,0,0,1]),1e-10); |
|
287 %! opts.droptol=1e-10; |
|
288 %! [l,u]=luinc(a,opts); |
|
289 %! assert(l*u, sparse([1,2,0,0;0,1,2,0;0,0,3,0;0,0,0,1]),1e-10); |
|
290 |
|
291 %!test |
|
292 %! a=sparse([1i,2,0,0;0,1,2,0;1e-14,0,3,0;0,0,0,1]); |
|
293 %! [l,u]=luinc(a,1e-10); |
|
294 %! assert(l*u, sparse([1i,2,0,0;0,1,2,0;0,0,3,0;0,0,0,1]),1e-10); |
|
295 %! opts.droptol=1e-10; |
|
296 %! [l,u]=luinc(a,opts); |
|
297 %! assert(l*u, sparse([1i,2,0,0;0,1,2,0;0,0,3,0;0,0,0,1]),1e-10); |
|
298 |
|
299 %!error splu(sparse([1i,2,0,0;0,1,2,0;1e-14,0,3,0;0,0,0,1])); |
|
300 |
|
301 */ |
|
302 |
|
303 /* |
5282
|
304 ;;; Local Variables: *** |
|
305 ;;; mode: C++ *** |
|
306 ;;; End: *** |
|
307 */ |