1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397
|
/*
Copyright (C) 2005-2015 David Bateman
This file is part of Octave.
Octave is free software; you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation; either version 3 of the License, or (at your
option) any later version.
Octave is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with Octave; see the file COPYING. If not, see
<http://www.gnu.org/licenses/>.
*/
#ifdef HAVE_CONFIG_H
#include <config.h>
#endif
#include "defun.h"
#include "error.h"
#include "gripes.h"
#include "oct-obj.h"
#include "utils.h"
#include "oct-map.h"
#include "MatrixType.h"
#include "SparseCmplxLU.h"
#include "SparsedbleLU.h"
#include "ov-re-sparse.h"
#include "ov-cx-sparse.h"
// FIXME: Deprecated in 4.0 and should be removed in 4.4.
DEFUN (__luinc__, args, nargout,
"-*- texinfo -*-\n\
@deftypefn {Built-in Function} {[@var{L}, @var{U}, @var{P}, @var{Q}] =} luinc (@var{A}, '0')\n\
@deftypefnx {Built-in Function} {[@var{L}, @var{U}, @var{P}, @var{Q}] =} luinc (@var{A}, @var{droptol})\n\
@deftypefnx {Built-in Function} {[@var{L}, @var{U}, @var{P}, @var{Q}] =} luinc (@var{A}, @var{opts})\n\
@cindex LU decomposition\n\
Produce the incomplete LU@tie{}factorization of the sparse matrix @var{A}.\n\
\n\
Two types of incomplete factorization are possible, and the type\n\
is determined by the second argument to @code{luinc}.\n\
\n\
Called with a second argument of @qcode{'0'}, the zero-level incomplete\n\
LU@tie{}factorization is produced. This creates a factorization of @var{A}\n\
where the position of the nonzero arguments correspond to the same\n\
positions as in the matrix @var{A}.\n\
\n\
Alternatively, the fill-in of the incomplete LU@tie{}factorization can\n\
be controlled through the variable @var{droptol} or the structure\n\
@var{opts}. The @sc{umfpack} multifrontal factorization code by Tim A.\n\
Davis is used for the incomplete LU@tie{}factorization, (availability\n\
@url{http://www.cise.ufl.edu/research/sparse/umfpack/})\n\
\n\
@var{droptol} determines the values below which the values in the\n\
LU@tie{} factorization are dropped and replaced by zero. It must be a\n\
positive scalar, and any values in the factorization whose absolute value\n\
are less than this value are dropped, expect if leaving them increase the\n\
sparsity of the matrix. Setting @var{droptol} to zero results in a complete\n\
LU@tie{}factorization which is the default.\n\
\n\
@var{opts} is a structure containing one or more of the fields\n\
\n\
@table @code\n\
@item droptol\n\
The drop tolerance as above. If @var{opts} only contains @code{droptol}\n\
then this is equivalent to using the variable @var{droptol}.\n\
\n\
@item milu\n\
A logical variable flagging whether to use the modified incomplete\n\
LU@tie{} factorization. In the case that @code{milu} is true, the dropped\n\
values are subtracted from the diagonal of the matrix @var{U} of the\n\
factorization. The default is @code{false}.\n\
\n\
@item udiag\n\
A logical variable that flags whether zero elements on the diagonal of\n\
@var{U} should be replaced with @var{droptol} to attempt to avoid singular\n\
factors. The default is @code{false}.\n\
\n\
@item thresh\n\
Defines the pivot threshold in the interval [0,1]. Values outside that\n\
range are ignored.\n\
@end table\n\
\n\
All other fields in @var{opts} are ignored. The outputs from @code{luinc}\n\
are the same as for @code{lu}.\n\
\n\
Given the string argument @qcode{\"vector\"}, @code{luinc} returns the\n\
values of @var{p} @var{q} as vector values.\n\
@seealso{sparse, lu, ilu, ichol}\n\
@end deftypefn")
{
int nargin = args.length ();
octave_value_list retval;
if (nargin == 0)
print_usage ();
else if (nargin < 2 || nargin > 3)
error ("luinc: incorrect number of arguments");
else
{
bool zero_level = false;
bool milu = false;
bool udiag = false;
Matrix thresh;
double droptol = -1.;
bool vecout = false;
if (args(1).is_string ())
{
if (args(1).string_value () == "0")
zero_level = true;
else
error ("luinc: unrecognized string argument");
}
else if (args(1).is_map ())
{
octave_scalar_map map = args(1).scalar_map_value ();
if (! error_state)
{
octave_value tmp;
tmp = map.getfield ("droptol");
if (tmp.is_defined ())
droptol = tmp.double_value ();
tmp = map.getfield ("milu");
if (tmp.is_defined ())
{
double val = tmp.double_value ();
milu = (val == 0. ? false : true);
}
tmp = map.getfield ("udiag");
if (tmp.is_defined ())
{
double val = tmp.double_value ();
udiag = (val == 0. ? false : true);
}
tmp = map.getfield ("thresh");
if (tmp.is_defined ())
{
thresh = tmp.matrix_value ();
if (thresh.nelem () == 1)
{
thresh.resize (1,2);
thresh(1) = thresh(0);
}
else if (thresh.nelem () != 2)
{
error ("luinc: expecting 2-element vector for thresh");
return retval;
}
}
}
else
{
error ("luinc: OPTS must be a scalar structure");
return retval;
}
}
else
droptol = args(1).double_value ();
if (nargin == 3)
{
std::string tmp = args(2).string_value ();
if (! error_state)
{
if (tmp.compare ("vector") == 0)
vecout = true;
else
error ("luinc: unrecognized string argument");
}
}
// FIXME: Add code for zero-level factorization
if (zero_level)
error ("luinc: zero-level factorization not implemented");
if (!error_state)
{
if (args(0).type_name () == "sparse matrix")
{
SparseMatrix sm = args(0).sparse_matrix_value ();
octave_idx_type sm_nr = sm.rows ();
octave_idx_type sm_nc = sm.cols ();
ColumnVector Qinit (sm_nc);
for (octave_idx_type i = 0; i < sm_nc; i++)
Qinit (i) = i;
if (! error_state)
{
switch (nargout)
{
case 0:
case 1:
case 2:
{
SparseLU fact (sm, Qinit, thresh, false, true, droptol,
milu, udiag);
if (! error_state)
{
SparseMatrix P = fact.Pr ();
SparseMatrix L = P.transpose () * fact.L ();
retval(1)
= octave_value (fact.U (),
MatrixType (MatrixType::Upper));
retval(0)
= octave_value (L, MatrixType
(MatrixType::Permuted_Lower,
sm_nr, fact.row_perm ()));
}
}
break;
case 3:
{
SparseLU fact (sm, Qinit, thresh, false, true, droptol,
milu, udiag);
if (! error_state)
{
if (vecout)
retval(2) = fact.Pr_vec ();
else
retval(2) = fact.Pr_mat ();
retval(1)
= octave_value (fact.U (),
MatrixType (MatrixType::Upper));
retval(0)
= octave_value (fact.L (),
MatrixType (MatrixType::Lower));
}
}
break;
case 4:
default:
{
SparseLU fact (sm, Qinit, thresh, false, false, droptol,
milu, udiag);
if (! error_state)
{
if (vecout)
{
retval(3) = fact.Pc_vec ();
retval(2) = fact.Pr_vec ();
}
else
{
retval(3) = fact.Pc_mat ();
retval(2) = fact.Pr_mat ();
}
retval(1)
= octave_value (fact.U (),
MatrixType (MatrixType::Upper));
retval(0)
= octave_value (fact.L (),
MatrixType (MatrixType::Lower));
}
}
break;
}
}
}
else if (args(0).type_name () == "sparse complex matrix")
{
SparseComplexMatrix sm =
args(0).sparse_complex_matrix_value ();
octave_idx_type sm_nr = sm.rows ();
octave_idx_type sm_nc = sm.cols ();
ColumnVector Qinit (sm_nc);
for (octave_idx_type i = 0; i < sm_nc; i++)
Qinit (i) = i;
if (! error_state)
{
switch (nargout)
{
case 0:
case 1:
case 2:
{
SparseComplexLU fact (sm, Qinit, thresh, false, true,
droptol, milu, udiag);
if (! error_state)
{
SparseMatrix P = fact.Pr ();
SparseComplexMatrix L = P.transpose () * fact.L ();
retval(1)
= octave_value (fact.U (),
MatrixType (MatrixType::Upper));
retval(0)
= octave_value (L, MatrixType
(MatrixType::Permuted_Lower,
sm_nr, fact.row_perm ()));
}
}
break;
case 3:
{
SparseComplexLU fact (sm, Qinit, thresh, false, true,
droptol, milu, udiag);
if (! error_state)
{
if (vecout)
retval(2) = fact.Pr_vec ();
else
retval(2) = fact.Pr_mat ();
retval(1)
= octave_value (fact.U (),
MatrixType (MatrixType::Upper));
retval(0)
= octave_value (fact.L (),
MatrixType (MatrixType::Lower));
}
}
break;
case 4:
default:
{
SparseComplexLU fact (sm, Qinit, thresh, false, false,
droptol, milu, udiag);
if (! error_state)
{
if (vecout)
{
retval(3) = fact.Pc_vec ();
retval(2) = fact.Pr_vec ();
}
else
{
retval(3) = fact.Pc_mat ();
retval(2) = fact.Pr_mat ();
}
retval(1)
= octave_value (fact.U (),
MatrixType (MatrixType::Upper));
retval(0)
= octave_value (fact.L (),
MatrixType (MatrixType::Lower));
}
}
break;
}
}
}
else
error ("luinc: matrix A must be sparse");
}
}
return retval;
}
/*
%!testif HAVE_UMFPACK
%! a = sparse ([1,2,0,0;0,1,2,0;1e-14,0,3,0;0,0,0,1]);
%! [l,u] = luinc (a, 1e-10);
%! assert (l*u, sparse ([1,2,0,0;0,1,2,0;0,0,3,0;0,0,0,1]), 1e-10);
%! opts.droptol = 1e-10;
%! [l,u] = luinc (a, opts);
%! assert (l*u, sparse ([1,2,0,0;0,1,2,0;0,0,3,0;0,0,0,1]), 1e-10);
%!testif HAVE_UMFPACK
%! a = sparse ([1i,2,0,0;0,1,2,0;1e-14,0,3,0;0,0,0,1]);
%! [l,u] = luinc (a, 1e-10);
%! assert (l*u, sparse ([1i,2,0,0;0,1,2,0;0,0,3,0;0,0,0,1]), 1e-10);
%! opts.droptol = 1e-10;
%! [l,u] = luinc (a, opts);
%! assert (l*u, sparse ([1i,2,0,0;0,1,2,0;0,0,3,0;0,0,0,1]), 1e-10);
*/
|