00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023 #ifdef HAVE_CONFIG_H
00024 #include <config.h>
00025 #endif
00026
00027 #include "defun-dld.h"
00028 #include "error.h"
00029 #include "gripes.h"
00030 #include "oct-obj.h"
00031 #include "utils.h"
00032 #include "oct-map.h"
00033
00034 #include "MatrixType.h"
00035 #include "SparseCmplxLU.h"
00036 #include "SparsedbleLU.h"
00037 #include "ov-re-sparse.h"
00038 #include "ov-cx-sparse.h"
00039
00040 DEFUN_DLD (luinc, args, nargout,
00041 "-*- texinfo -*-\n\
00042 @deftypefn {Loadable Function} {[@var{L}, @var{U}, @var{P}, @var{Q}] =} luinc (@var{A}, '0')\n\
00043 @deftypefnx {Loadable Function} {[@var{L}, @var{U}, @var{P}, @var{Q}] =} luinc (@var{A}, @var{droptol})\n\
00044 @deftypefnx {Loadable Function} {[@var{L}, @var{U}, @var{P}, @var{Q}] =} luinc (@var{A}, @var{opts})\n\
00045 @cindex LU decomposition\n\
00046 Produce the incomplete LU@tie{}factorization of the sparse matrix @var{A}.\n\
00047 Two types of incomplete factorization are possible, and the type\n\
00048 is determined by the second argument to @code{luinc}.\n\
00049 \n\
00050 Called with a second argument of '0', the zero-level incomplete\n\
00051 LU@tie{}factorization is produced. This creates a factorization of @var{A}\n\
00052 where the position of the non-zero arguments correspond to the same\n\
00053 positions as in the matrix @var{A}.\n\
00054 \n\
00055 Alternatively, the fill-in of the incomplete LU@tie{}factorization can\n\
00056 be controlled through the variable @var{droptol} or the structure\n\
00057 @var{opts}. The @sc{umfpack} multifrontal factorization code by Tim A.\n\
00058 Davis is used for the incomplete LU@tie{}factorization, (availability\n\
00059 @url{http://www.cise.ufl.edu/research/sparse/umfpack/})\n\
00060 \n\
00061 @var{droptol} determines the values below which the values in the\n\
00062 LU@tie{} factorization are dropped and replaced by zero. It must be a\n\
00063 positive scalar, and any values in the factorization whose absolute value\n\
00064 are less than this value are dropped, expect if leaving them increase the\n\
00065 sparsity of the matrix. Setting @var{droptol} to zero results in a complete\n\
00066 LU@tie{}factorization which is the default.\n\
00067 \n\
00068 @var{opts} is a structure containing one or more of the fields\n\
00069 \n\
00070 @table @code\n\
00071 @item droptol\n\
00072 The drop tolerance as above. If @var{opts} only contains @code{droptol}\n\
00073 then this is equivalent to using the variable @var{droptol}.\n\
00074 \n\
00075 @item milu\n\
00076 A logical variable flagging whether to use the modified incomplete\n\
00077 LU@tie{} factorization. In the case that @code{milu} is true, the dropped\n\
00078 values are subtracted from the diagonal of the matrix @var{U} of the\n\
00079 factorization. The default is @code{false}.\n\
00080 \n\
00081 @item udiag\n\
00082 A logical variable that flags whether zero elements on the diagonal of\n\
00083 @var{U} should be replaced with @var{droptol} to attempt to avoid singular\n\
00084 factors. The default is @code{false}.\n\
00085 \n\
00086 @item thresh\n\
00087 Defines the pivot threshold in the interval [0,1]. Values outside that\n\
00088 range are ignored.\n\
00089 @end table\n\
00090 \n\
00091 All other fields in @var{opts} are ignored. The outputs from @code{luinc}\n\
00092 are the same as for @code{lu}.\n\
00093 \n\
00094 Given the string argument 'vector', @code{luinc} returns the values of\n\
00095 @var{p} @var{q} as vector values.\n\
00096 @seealso{sparse, lu}\n\
00097 @end deftypefn")
00098 {
00099 int nargin = args.length ();
00100 octave_value_list retval;
00101
00102 if (nargin == 0)
00103 print_usage ();
00104 else if (nargin < 2 || nargin > 3)
00105 error ("luinc: incorrect number of arguments");
00106 else
00107 {
00108 bool zero_level = false;
00109 bool milu = false;
00110 bool udiag = false;
00111 Matrix thresh;
00112 double droptol = -1.;
00113 bool vecout = false;
00114
00115 if (args(1).is_string ())
00116 {
00117 if (args(1).string_value () == "0")
00118 zero_level = true;
00119 else
00120 error ("luinc: unrecognized string argument");
00121 }
00122 else if (args(1).is_map ())
00123 {
00124 octave_scalar_map map = args(1).scalar_map_value ();
00125
00126 if (! error_state)
00127 {
00128 octave_value tmp;
00129
00130 tmp = map.getfield ("droptol");
00131 if (tmp.is_defined ())
00132 droptol = tmp.double_value ();
00133
00134 tmp = map.getfield ("milu");
00135 if (tmp.is_defined ())
00136 {
00137 double val = tmp.double_value ();
00138
00139 milu = (val == 0. ? false : true);
00140 }
00141
00142 tmp = map.getfield ("udiag");
00143 if (tmp.is_defined ())
00144 {
00145 double val = tmp.double_value ();
00146
00147 udiag = (val == 0. ? false : true);
00148 }
00149
00150 tmp = map.getfield ("thresh");
00151 if (tmp.is_defined ())
00152 {
00153 thresh = tmp.matrix_value ();
00154
00155 if (thresh.nelem () == 1)
00156 {
00157 thresh.resize(1,2);
00158 thresh(1) = thresh(0);
00159 }
00160 else if (thresh.nelem () != 2)
00161 {
00162 error ("luinc: expecting 2-element vector for thresh");
00163 return retval;
00164 }
00165 }
00166 }
00167 else
00168 {
00169 error ("luinc: OPTS must be a scalar structure");
00170 return retval;
00171 }
00172 }
00173 else
00174 droptol = args(1).double_value ();
00175
00176 if (nargin == 3)
00177 {
00178 std::string tmp = args(2).string_value ();
00179
00180 if (! error_state )
00181 {
00182 if (tmp.compare ("vector") == 0)
00183 vecout = true;
00184 else
00185 error ("luinc: unrecognized string argument");
00186 }
00187 }
00188
00189
00190 if (zero_level)
00191 error ("luinc: zero-level factorization not implemented");
00192
00193 if (!error_state)
00194 {
00195 if (args(0).type_name () == "sparse matrix")
00196 {
00197 SparseMatrix sm = args(0).sparse_matrix_value ();
00198 octave_idx_type sm_nr = sm.rows ();
00199 octave_idx_type sm_nc = sm.cols ();
00200 ColumnVector Qinit (sm_nc);
00201
00202 for (octave_idx_type i = 0; i < sm_nc; i++)
00203 Qinit (i) = i;
00204
00205 if (! error_state)
00206 {
00207 switch (nargout)
00208 {
00209 case 0:
00210 case 1:
00211 case 2:
00212 {
00213 SparseLU fact (sm, Qinit, thresh, false, true, droptol,
00214 milu, udiag);
00215
00216 if (! error_state)
00217 {
00218 SparseMatrix P = fact.Pr ();
00219 SparseMatrix L = P.transpose () * fact.L ();
00220 retval(1) = octave_value (fact.U (),
00221 MatrixType (MatrixType::Upper));
00222 retval(0) = octave_value (L, MatrixType
00223 (MatrixType::Permuted_Lower,
00224 sm_nr, fact.row_perm ()));
00225 }
00226 }
00227 break;
00228
00229 case 3:
00230 {
00231 SparseLU fact (sm, Qinit, thresh, false, true, droptol,
00232 milu, udiag);
00233
00234 if (! error_state)
00235 {
00236 if (vecout)
00237 retval(2) = fact.Pr_vec ();
00238 else
00239 retval(2) = fact.Pr_mat ();
00240 retval(1) = octave_value (fact.U (),
00241 MatrixType (MatrixType::Upper));
00242 retval(0) = octave_value (fact.L (),
00243 MatrixType (MatrixType::Lower));
00244 }
00245 }
00246 break;
00247
00248 case 4:
00249 default:
00250 {
00251 SparseLU fact (sm, Qinit, thresh, false, false, droptol,
00252 milu, udiag);
00253
00254 if (! error_state)
00255 {
00256 if (vecout)
00257 {
00258 retval(3) = fact.Pc_vec ();
00259 retval(2) = fact.Pr_vec ();
00260 }
00261 else
00262 {
00263 retval(3) = fact.Pc_mat ();
00264 retval(2) = fact.Pr_mat ();
00265 }
00266 retval(1) = octave_value (fact.U (),
00267 MatrixType (MatrixType::Upper));
00268 retval(0) = octave_value (fact.L (),
00269 MatrixType (MatrixType::Lower));
00270 }
00271 }
00272 break;
00273 }
00274 }
00275 }
00276 else if (args(0).type_name () == "sparse complex matrix")
00277 {
00278 SparseComplexMatrix sm =
00279 args(0).sparse_complex_matrix_value ();
00280 octave_idx_type sm_nr = sm.rows ();
00281 octave_idx_type sm_nc = sm.cols ();
00282 ColumnVector Qinit (sm_nc);
00283
00284 for (octave_idx_type i = 0; i < sm_nc; i++)
00285 Qinit (i) = i;
00286
00287 if (! error_state)
00288 {
00289 switch (nargout)
00290 {
00291 case 0:
00292 case 1:
00293 case 2:
00294 {
00295 SparseComplexLU fact (sm, Qinit, thresh, false, true,
00296 droptol, milu, udiag);
00297
00298
00299 if (! error_state)
00300 {
00301 SparseMatrix P = fact.Pr ();
00302 SparseComplexMatrix L = P.transpose () * fact.L ();
00303 retval(1) = octave_value (fact.U (),
00304 MatrixType (MatrixType::Upper));
00305 retval(0) = octave_value (L, MatrixType
00306 (MatrixType::Permuted_Lower,
00307 sm_nr, fact.row_perm ()));
00308 }
00309 }
00310 break;
00311
00312 case 3:
00313 {
00314 SparseComplexLU fact (sm, Qinit, thresh, false, true,
00315 droptol, milu, udiag);
00316
00317 if (! error_state)
00318 {
00319 if (vecout)
00320 retval(2) = fact.Pr_vec ();
00321 else
00322 retval(2) = fact.Pr_mat ();
00323 retval(1) = octave_value (fact.U (),
00324 MatrixType (MatrixType::Upper));
00325 retval(0) = octave_value (fact.L (),
00326 MatrixType (MatrixType::Lower));
00327 }
00328 }
00329 break;
00330
00331 case 4:
00332 default:
00333 {
00334 SparseComplexLU fact (sm, Qinit, thresh, false, false,
00335 droptol, milu, udiag);
00336
00337 if (! error_state)
00338 {
00339 if (vecout)
00340 {
00341 retval(3) = fact.Pc_vec ();
00342 retval(2) = fact.Pr_vec ();
00343 }
00344 else
00345 {
00346 retval(3) = fact.Pc_mat ();
00347 retval(2) = fact.Pr_mat ();
00348 }
00349 retval(1) = octave_value (fact.U (),
00350 MatrixType (MatrixType::Upper));
00351 retval(0) = octave_value (fact.L (),
00352 MatrixType (MatrixType::Lower));
00353 }
00354 }
00355 break;
00356 }
00357 }
00358 }
00359 else
00360 error ("luinc: matrix A must be sparse");
00361 }
00362 }
00363
00364 return retval;
00365 }
00366
00367
00368
00369
00370
00371
00372
00373
00374
00375
00376
00377
00378
00379
00380
00381
00382
00383
00384
00385