GNU Octave  3.8.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
luinc.cc
Go to the documentation of this file.
1 /*
2 
3 Copyright (C) 2005-2013 David Bateman
4 
5 This file is part of Octave.
6 
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
9 Free Software Foundation; either version 3 of the License, or (at your
10 option) any later version.
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
18 along with Octave; see the file COPYING. If not, see
19 <http://www.gnu.org/licenses/>.
20 
21 */
22 
23 #ifdef HAVE_CONFIG_H
24 #include <config.h>
25 #endif
26 
27 #include "defun.h"
28 #include "error.h"
29 #include "gripes.h"
30 #include "oct-obj.h"
31 #include "utils.h"
32 #include "oct-map.h"
33 
34 #include "MatrixType.h"
35 #include "SparseCmplxLU.h"
36 #include "SparsedbleLU.h"
37 #include "ov-re-sparse.h"
38 #include "ov-cx-sparse.h"
39 
40 DEFUN (luinc, args, nargout,
41  "-*- texinfo -*-\n\
42 @deftypefn {Built-in Function} {[@var{L}, @var{U}, @var{P}, @var{Q}] =} luinc (@var{A}, '0')\n\
43 @deftypefnx {Built-in Function} {[@var{L}, @var{U}, @var{P}, @var{Q}] =} luinc (@var{A}, @var{droptol})\n\
44 @deftypefnx {Built-in 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@tie{}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 @code{luinc}.\n\
49 \n\
50 Called with a second argument of @qcode{'0'}, the zero-level incomplete\n\
51 LU@tie{}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@tie{}factorization can\n\
56 be controlled through the variable @var{droptol} or the structure\n\
57 @var{opts}. The @sc{umfpack} multifrontal factorization code by Tim A.\n\
58 Davis is used for the incomplete LU@tie{}factorization, (availability\n\
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\n\
62 LU@tie{} factorization are dropped and replaced by zero. It must be a\n\
63 positive scalar, and any values in the factorization whose absolute value\n\
64 are 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 complete\n\
66 LU@tie{}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\n\
77 LU@tie{} factorization. In the case that @code{milu} is true, the dropped\n\
78 values are subtracted from the diagonal of the matrix @var{U} of the\n\
79 factorization. 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\n\
83 @var{U} 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 @code{luinc}\n\
92 are the same as for @code{lu}.\n\
93 \n\
94 Given the string argument @qcode{\"vector\"}, @code{luinc} returns the\n\
95 values of @var{p} @var{q} as vector values.\n\
96 @seealso{sparse, lu}\n\
97 @end deftypefn")
98 {
99  int nargin = args.length ();
100  octave_value_list retval;
101 
102  if (nargin == 0)
103  print_usage ();
104  else if (nargin < 2 || nargin > 3)
105  error ("luinc: incorrect number of arguments");
106  else
107  {
108  bool zero_level = false;
109  bool milu = false;
110  bool udiag = false;
111  Matrix thresh;
112  double droptol = -1.;
113  bool vecout = false;
114 
115  if (args(1).is_string ())
116  {
117  if (args(1).string_value () == "0")
118  zero_level = true;
119  else
120  error ("luinc: unrecognized string argument");
121  }
122  else if (args(1).is_map ())
123  {
124  octave_scalar_map map = args(1).scalar_map_value ();
125 
126  if (! error_state)
127  {
128  octave_value tmp;
129 
130  tmp = map.getfield ("droptol");
131  if (tmp.is_defined ())
132  droptol = tmp.double_value ();
133 
134  tmp = map.getfield ("milu");
135  if (tmp.is_defined ())
136  {
137  double val = tmp.double_value ();
138 
139  milu = (val == 0. ? false : true);
140  }
141 
142  tmp = map.getfield ("udiag");
143  if (tmp.is_defined ())
144  {
145  double val = tmp.double_value ();
146 
147  udiag = (val == 0. ? false : true);
148  }
149 
150  tmp = map.getfield ("thresh");
151  if (tmp.is_defined ())
152  {
153  thresh = tmp.matrix_value ();
154 
155  if (thresh.nelem () == 1)
156  {
157  thresh.resize (1,2);
158  thresh(1) = thresh(0);
159  }
160  else if (thresh.nelem () != 2)
161  {
162  error ("luinc: expecting 2-element vector for thresh");
163  return retval;
164  }
165  }
166  }
167  else
168  {
169  error ("luinc: OPTS must be a scalar structure");
170  return retval;
171  }
172  }
173  else
174  droptol = args(1).double_value ();
175 
176  if (nargin == 3)
177  {
178  std::string tmp = args(2).string_value ();
179 
180  if (! error_state )
181  {
182  if (tmp.compare ("vector") == 0)
183  vecout = true;
184  else
185  error ("luinc: unrecognized string argument");
186  }
187  }
188 
189  // FIXME: Add code for zero-level factorization
190  if (zero_level)
191  error ("luinc: zero-level factorization not implemented");
192 
193  if (!error_state)
194  {
195  if (args(0).type_name () == "sparse matrix")
196  {
197  SparseMatrix sm = args(0).sparse_matrix_value ();
198  octave_idx_type sm_nr = sm.rows ();
199  octave_idx_type sm_nc = sm.cols ();
200  ColumnVector Qinit (sm_nc);
201 
202  for (octave_idx_type i = 0; i < sm_nc; i++)
203  Qinit (i) = i;
204 
205  if (! error_state)
206  {
207  switch (nargout)
208  {
209  case 0:
210  case 1:
211  case 2:
212  {
213  SparseLU fact (sm, Qinit, thresh, false, true, droptol,
214  milu, udiag);
215 
216  if (! error_state)
217  {
218  SparseMatrix P = fact.Pr ();
219  SparseMatrix L = P.transpose () * fact.L ();
220  retval(1)
221  = octave_value (fact.U (),
223  retval(0)
226  sm_nr, fact.row_perm ()));
227  }
228  }
229  break;
230 
231  case 3:
232  {
233  SparseLU fact (sm, Qinit, thresh, false, true, droptol,
234  milu, udiag);
235 
236  if (! error_state)
237  {
238  if (vecout)
239  retval(2) = fact.Pr_vec ();
240  else
241  retval(2) = fact.Pr_mat ();
242  retval(1)
243  = octave_value (fact.U (),
245  retval(0)
246  = octave_value (fact.L (),
248  }
249  }
250  break;
251 
252  case 4:
253  default:
254  {
255  SparseLU fact (sm, Qinit, thresh, false, false, droptol,
256  milu, udiag);
257 
258  if (! error_state)
259  {
260  if (vecout)
261  {
262  retval(3) = fact.Pc_vec ();
263  retval(2) = fact.Pr_vec ();
264  }
265  else
266  {
267  retval(3) = fact.Pc_mat ();
268  retval(2) = fact.Pr_mat ();
269  }
270  retval(1)
271  = octave_value (fact.U (),
273  retval(0)
274  = octave_value (fact.L (),
276  }
277  }
278  break;
279  }
280  }
281  }
282  else if (args(0).type_name () == "sparse complex matrix")
283  {
285  args(0).sparse_complex_matrix_value ();
286  octave_idx_type sm_nr = sm.rows ();
287  octave_idx_type sm_nc = sm.cols ();
288  ColumnVector Qinit (sm_nc);
289 
290  for (octave_idx_type i = 0; i < sm_nc; i++)
291  Qinit (i) = i;
292 
293  if (! error_state)
294  {
295  switch (nargout)
296  {
297  case 0:
298  case 1:
299  case 2:
300  {
301  SparseComplexLU fact (sm, Qinit, thresh, false, true,
302  droptol, milu, udiag);
303 
304 
305  if (! error_state)
306  {
307  SparseMatrix P = fact.Pr ();
308  SparseComplexMatrix L = P.transpose () * fact.L ();
309  retval(1)
310  = octave_value (fact.U (),
312  retval(0)
315  sm_nr, fact.row_perm ()));
316  }
317  }
318  break;
319 
320  case 3:
321  {
322  SparseComplexLU fact (sm, Qinit, thresh, false, true,
323  droptol, milu, udiag);
324 
325  if (! error_state)
326  {
327  if (vecout)
328  retval(2) = fact.Pr_vec ();
329  else
330  retval(2) = fact.Pr_mat ();
331  retval(1)
332  = octave_value (fact.U (),
334  retval(0)
335  = octave_value (fact.L (),
337  }
338  }
339  break;
340 
341  case 4:
342  default:
343  {
344  SparseComplexLU fact (sm, Qinit, thresh, false, false,
345  droptol, milu, udiag);
346 
347  if (! error_state)
348  {
349  if (vecout)
350  {
351  retval(3) = fact.Pc_vec ();
352  retval(2) = fact.Pr_vec ();
353  }
354  else
355  {
356  retval(3) = fact.Pc_mat ();
357  retval(2) = fact.Pr_mat ();
358  }
359  retval(1)
360  = octave_value (fact.U (),
362  retval(0)
363  = octave_value (fact.L (),
365  }
366  }
367  break;
368  }
369  }
370  }
371  else
372  error ("luinc: matrix A must be sparse");
373  }
374  }
375 
376  return retval;
377 }
378 
379 /*
380 %!testif HAVE_UMFPACK
381 %! a = sparse ([1,2,0,0;0,1,2,0;1e-14,0,3,0;0,0,0,1]);
382 %! [l,u] = luinc (a, 1e-10);
383 %! assert (l*u, sparse ([1,2,0,0;0,1,2,0;0,0,3,0;0,0,0,1]), 1e-10);
384 %! opts.droptol = 1e-10;
385 %! [l,u] = luinc (a, opts);
386 %! assert (l*u, sparse ([1,2,0,0;0,1,2,0;0,0,3,0;0,0,0,1]), 1e-10);
387 
388 %!testif HAVE_UMFPACK
389 %! a = sparse ([1i,2,0,0;0,1,2,0;1e-14,0,3,0;0,0,0,1]);
390 %! [l,u] = luinc (a, 1e-10);
391 %! assert (l*u, sparse ([1i,2,0,0;0,1,2,0;0,0,3,0;0,0,0,1]), 1e-10);
392 %! opts.droptol = 1e-10;
393 %! [l,u] = luinc (a, opts);
394 %! assert (l*u, sparse ([1i,2,0,0;0,1,2,0;0,0,3,0;0,0,0,1]), 1e-10);
395 */