Mercurial > octave
annotate libinterp/corefcn/lu.cc @ 31607:aac27ad79be6 stable
maint: Re-indent code after switch to using namespace macros.
* build-env.h, build-env.in.cc, Cell.h, __betainc__.cc, __eigs__.cc,
__ftp__.cc, __ichol__.cc, __ilu__.cc, __isprimelarge__.cc, __magick_read__.cc,
__pchip_deriv__.cc, amd.cc, base-text-renderer.cc, base-text-renderer.h,
besselj.cc, bitfcns.cc, bsxfun.cc, c-file-ptr-stream.h, call-stack.cc,
call-stack.h, ccolamd.cc, cellfun.cc, chol.cc, colamd.cc, dasrt.cc, data.cc,
debug.cc, defaults.cc, defaults.h, det.cc, display.cc, display.h, dlmread.cc,
dynamic-ld.cc, dynamic-ld.h, ellipj.cc, environment.cc, environment.h,
error.cc, error.h, errwarn.h, event-manager.cc, event-manager.h,
event-queue.cc, event-queue.h, fcn-info.cc, fcn-info.h, fft.cc, fft2.cc,
file-io.cc, filter.cc, find.cc, ft-text-renderer.cc, ft-text-renderer.h,
gcd.cc, gl-render.cc, gl-render.h, gl2ps-print.cc, gl2ps-print.h,
graphics-toolkit.cc, graphics-toolkit.h, graphics.cc, gsvd.cc, gtk-manager.cc,
gtk-manager.h, help.cc, help.h, hook-fcn.cc, hook-fcn.h, input.cc, input.h,
interpreter-private.cc, interpreter-private.h, interpreter.cc, interpreter.h,
inv.cc, jsondecode.cc, jsonencode.cc, latex-text-renderer.cc,
latex-text-renderer.h, load-path.cc, load-path.h, load-save.cc, load-save.h,
lookup.cc, ls-hdf5.cc, ls-mat4.cc, ls-mat5.cc, lsode.cc, lu.cc, mappers.cc,
matrix_type.cc, max.cc, mex.cc, mexproto.h, mxarray.h, mxtypes.in.h,
oct-errno.in.cc, oct-hdf5-types.cc, oct-hist.cc, oct-hist.h, oct-map.cc,
oct-map.h, oct-opengl.h, oct-prcstrm.h, oct-process.cc, oct-process.h,
oct-stdstrm.h, oct-stream.cc, oct-stream.h, oct-strstrm.h,
octave-default-image.h, ordqz.cc, ordschur.cc, pager.cc, pager.h, pinv.cc,
pow2.cc, pr-output.cc, psi.cc, qr.cc, quadcc.cc, rand.cc, regexp.cc,
settings.cc, settings.h, sighandlers.cc, sighandlers.h, sparse-xpow.cc,
sqrtm.cc, stack-frame.cc, stack-frame.h, stream-euler.cc, strfns.cc, svd.cc,
syminfo.cc, syminfo.h, symrcm.cc, symrec.cc, symrec.h, symscope.cc, symscope.h,
symtab.cc, symtab.h, sysdep.cc, sysdep.h, text-engine.cc, text-engine.h,
text-renderer.cc, text-renderer.h, time.cc, toplev.cc, typecast.cc,
url-handle-manager.cc, url-handle-manager.h, urlwrite.cc, utils.cc, utils.h,
variables.cc, variables.h, xdiv.cc, __delaunayn__.cc, __init_fltk__.cc,
__init_gnuplot__.cc, __ode15__.cc, __voronoi__.cc, audioread.cc, convhulln.cc,
gzip.cc, cdef-class.cc, cdef-class.h, cdef-fwd.h, cdef-manager.cc,
cdef-manager.h, cdef-method.cc, cdef-method.h, cdef-object.cc, cdef-object.h,
cdef-package.cc, cdef-package.h, cdef-property.cc, cdef-property.h,
cdef-utils.cc, cdef-utils.h, ov-base-diag.cc, ov-base-int.cc, ov-base-mat.cc,
ov-base-mat.h, ov-base-scalar.cc, ov-base.cc, ov-base.h, ov-bool-mat.cc,
ov-bool-mat.h, ov-bool-sparse.cc, ov-bool.cc, ov-builtin.h, ov-cell.cc,
ov-ch-mat.cc, ov-class.cc, ov-class.h, ov-classdef.cc, ov-classdef.h,
ov-complex.cc, ov-cx-diag.cc, ov-cx-mat.cc, ov-cx-sparse.cc, ov-dld-fcn.cc,
ov-dld-fcn.h, ov-fcn-handle.cc, ov-fcn-handle.h, ov-fcn.h, ov-float.cc,
ov-flt-complex.cc, ov-flt-cx-diag.cc, ov-flt-cx-mat.cc, ov-flt-re-diag.cc,
ov-flt-re-mat.cc, ov-flt-re-mat.h, ov-intx.h, ov-java.cc, ov-lazy-idx.cc,
ov-legacy-range.cc, ov-magic-int.cc, ov-mex-fcn.cc, ov-mex-fcn.h,
ov-null-mat.cc, ov-perm.cc, ov-range.cc, ov-re-diag.cc, ov-re-mat.cc,
ov-re-mat.h, ov-re-sparse.cc, ov-scalar.cc, ov-str-mat.cc, ov-struct.cc,
ov-typeinfo.cc, ov-typeinfo.h, ov-usr-fcn.cc, ov-usr-fcn.h, ov.cc, ov.h, ovl.h,
octave.cc, octave.h, op-b-sbm.cc, op-bm-sbm.cc, op-cs-scm.cc, op-fm-fcm.cc,
op-fs-fcm.cc, op-s-scm.cc, op-scm-cs.cc, op-scm-s.cc, op-sm-cs.cc, ops.h,
anon-fcn-validator.cc, anon-fcn-validator.h, bp-table.cc, bp-table.h,
comment-list.cc, comment-list.h, filepos.h, lex.h, oct-lvalue.cc, oct-lvalue.h,
parse.h, profiler.cc, profiler.h, pt-anon-scopes.cc, pt-anon-scopes.h,
pt-arg-list.cc, pt-arg-list.h, pt-args-block.cc, pt-args-block.h,
pt-array-list.cc, pt-array-list.h, pt-assign.cc, pt-assign.h, pt-binop.cc,
pt-binop.h, pt-bp.cc, pt-bp.h, pt-cbinop.cc, pt-cbinop.h, pt-cell.cc,
pt-cell.h, pt-check.cc, pt-check.h, pt-classdef.cc, pt-classdef.h, pt-cmd.h,
pt-colon.cc, pt-colon.h, pt-const.cc, pt-const.h, pt-decl.cc, pt-decl.h,
pt-eval.cc, pt-eval.h, pt-except.cc, pt-except.h, pt-exp.cc, pt-exp.h,
pt-fcn-handle.cc, pt-fcn-handle.h, pt-id.cc, pt-id.h, pt-idx.cc, pt-idx.h,
pt-jump.h, pt-loop.cc, pt-loop.h, pt-mat.cc, pt-mat.h, pt-misc.cc, pt-misc.h,
pt-pr-code.cc, pt-pr-code.h, pt-select.cc, pt-select.h, pt-spmd.cc, pt-spmd.h,
pt-stmt.cc, pt-stmt.h, pt-tm-const.cc, pt-tm-const.h, pt-unop.cc, pt-unop.h,
pt-walk.cc, pt-walk.h, pt.cc, pt.h, token.cc, token.h, Range.cc, Range.h,
idx-vector.cc, idx-vector.h, range-fwd.h, CollocWt.cc, CollocWt.h,
aepbalance.cc, aepbalance.h, chol.cc, chol.h, gepbalance.cc, gepbalance.h,
gsvd.cc, gsvd.h, hess.cc, hess.h, lo-mappers.cc, lo-mappers.h, lo-specfun.cc,
lo-specfun.h, lu.cc, lu.h, oct-convn.cc, oct-convn.h, oct-fftw.cc, oct-fftw.h,
oct-norm.cc, oct-norm.h, oct-rand.cc, oct-rand.h, oct-spparms.cc,
oct-spparms.h, qr.cc, qr.h, qrp.cc, qrp.h, randgamma.cc, randgamma.h,
randmtzig.cc, randmtzig.h, randpoisson.cc, randpoisson.h, schur.cc, schur.h,
sparse-chol.cc, sparse-chol.h, sparse-lu.cc, sparse-lu.h, sparse-qr.cc,
sparse-qr.h, svd.cc, svd.h, child-list.cc, child-list.h, dir-ops.cc, dir-ops.h,
file-ops.cc, file-ops.h, file-stat.cc, file-stat.h, lo-sysdep.cc, lo-sysdep.h,
lo-sysinfo.cc, lo-sysinfo.h, mach-info.cc, mach-info.h, oct-env.cc, oct-env.h,
oct-group.cc, oct-group.h, oct-password.cc, oct-password.h, oct-syscalls.cc,
oct-syscalls.h, oct-time.cc, oct-time.h, oct-uname.cc, oct-uname.h,
action-container.cc, action-container.h, base-list.h, cmd-edit.cc, cmd-edit.h,
cmd-hist.cc, cmd-hist.h, f77-fcn.h, file-info.cc, file-info.h,
lo-array-errwarn.cc, lo-array-errwarn.h, lo-hash.cc, lo-hash.h, lo-ieee.h,
lo-regexp.cc, lo-regexp.h, lo-utils.cc, lo-utils.h, oct-base64.cc,
oct-base64.h, oct-glob.cc, oct-glob.h, oct-inttypes.h, oct-mutex.cc,
oct-mutex.h, oct-refcount.h, oct-shlib.cc, oct-shlib.h, oct-sparse.cc,
oct-sparse.h, oct-string.h, octave-preserve-stream-state.h, pathsearch.cc,
pathsearch.h, quit.cc, quit.h, unwind-prot.cc, unwind-prot.h, url-transfer.cc,
url-transfer.h:
Re-indent code after switch to using namespace macros.
author | Rik <rik@octave.org> |
---|---|
date | Thu, 01 Dec 2022 18:02:15 -0800 |
parents | e88a07dec498 |
children | 597f3ee61a48 |
rev | line source |
---|---|
27923
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
1 //////////////////////////////////////////////////////////////////////// |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
2 // |
30564
796f54d4ddbf
update Octave Project Developers copyright for the new year
John W. Eaton <jwe@octave.org>
parents:
30390
diff
changeset
|
3 // Copyright (C) 1996-2022 The Octave Project Developers |
27923
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
4 // |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
5 // See the file COPYRIGHT.md in the top-level directory of this |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
6 // distribution or <https://octave.org/copyright/>. |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
7 // |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
8 // This file is part of Octave. |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
9 // |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
10 // Octave is free software: you can redistribute it and/or modify it |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
11 // under the terms of the GNU General Public License as published by |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
12 // the Free Software Foundation, either version 3 of the License, or |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
13 // (at your option) any later version. |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
14 // |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
15 // Octave is distributed in the hope that it will be useful, but |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
16 // WITHOUT ANY WARRANTY; without even the implied warranty of |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
17 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
18 // GNU General Public License for more details. |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
19 // |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
20 // You should have received a copy of the GNU General Public License |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
21 // along with Octave; see the file COPYING. If not, see |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
22 // <https://www.gnu.org/licenses/>. |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
23 // |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
24 //////////////////////////////////////////////////////////////////////// |
2928 | 25 |
21724
aba2e6293dd8
use "#if ..." consistently instead of "#ifdef" and "#ifndef"
John W. Eaton <jwe@octave.org>
parents:
21547
diff
changeset
|
26 #if defined (HAVE_CONFIG_H) |
21301
40de9f8f23a6
Use '#include "config.h"' rather than <config.h>.
Rik <rik@octave.org>
parents:
21271
diff
changeset
|
27 # include "config.h" |
2928 | 28 #endif |
29 | |
21271
7e67c7f82fc1
better use of templates for lu factorization classes
John W. Eaton <jwe@octave.org>
parents:
21200
diff
changeset
|
30 #include "lu.h" |
21146
ea9c05014809
revamp sparse LU factorization classes
John W. Eaton <jwe@octave.org>
parents:
21139
diff
changeset
|
31 #include "sparse-lu.h" |
2928 | 32 |
15039
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14854
diff
changeset
|
33 #include "defun.h" |
2928 | 34 #include "error.h" |
21100
e39e05d90788
Switch gripe_XXX to either err_XXX or warn_XXX naming scheme.
Rik <rik@octave.org>
parents:
21055
diff
changeset
|
35 #include "errwarn.h" |
20940
48b2ad5ee801
maint: Rename oct-obj.[cc|h] to ovl.[cc|h] for clarity.
Rik <rik@octave.org>
parents:
20921
diff
changeset
|
36 #include "ovl.h" |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
37 #include "ov-re-sparse.h" |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
38 #include "ov-cx-sparse.h" |
2928 | 39 |
31605
e88a07dec498
maint: Use macros to begin/end C++ namespaces.
Rik <rik@octave.org>
parents:
30565
diff
changeset
|
40 OCTAVE_BEGIN_NAMESPACE(octave) |
29958
32c3a5805893
move DEFUN and DEFMETHOD functions inside octave namespace
John W. Eaton <jwe@octave.org>
parents:
29359
diff
changeset
|
41 |
21139
538b57866b90
consistently use "typename" intead of "class" in template declarations
John W. Eaton <jwe@octave.org>
parents:
21100
diff
changeset
|
42 template <typename MT> |
8869
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
43 static octave_value |
29961
7d6709900da7
eliminate octave:: namespace tags in DEFUN and DEFMETHOD and more
John W. Eaton <jwe@octave.org>
parents:
29958
diff
changeset
|
44 get_lu_l (const math::lu<MT>& fact) |
8869
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
45 { |
9715
9f27172fbd1e
auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents:
9709
diff
changeset
|
46 MT L = fact.L (); |
23593
a8361bc2361a
maint: Deprecate is_square and replace with issquare.
Rik <rik@octave.org>
parents:
23586
diff
changeset
|
47 if (L.issquare ()) |
9715
9f27172fbd1e
auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents:
9709
diff
changeset
|
48 return octave_value (L, MatrixType (MatrixType::Lower)); |
8869
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
49 else |
9715
9f27172fbd1e
auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents:
9709
diff
changeset
|
50 return L; |
9f27172fbd1e
auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents:
9709
diff
changeset
|
51 } |
8869
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
52 |
21139
538b57866b90
consistently use "typename" intead of "class" in template declarations
John W. Eaton <jwe@octave.org>
parents:
21100
diff
changeset
|
53 template <typename MT> |
9715
9f27172fbd1e
auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents:
9709
diff
changeset
|
54 static octave_value |
29961
7d6709900da7
eliminate octave:: namespace tags in DEFUN and DEFMETHOD and more
John W. Eaton <jwe@octave.org>
parents:
29958
diff
changeset
|
55 get_lu_u (const math::lu<MT>& fact) |
9715
9f27172fbd1e
auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents:
9709
diff
changeset
|
56 { |
9f27172fbd1e
auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents:
9709
diff
changeset
|
57 MT U = fact.U (); |
23593
a8361bc2361a
maint: Deprecate is_square and replace with issquare.
Rik <rik@octave.org>
parents:
23586
diff
changeset
|
58 if (U.issquare () && fact.regular ()) |
9715
9f27172fbd1e
auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents:
9709
diff
changeset
|
59 return octave_value (U, MatrixType (MatrixType::Upper)); |
9f27172fbd1e
auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents:
9709
diff
changeset
|
60 else |
9f27172fbd1e
auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents:
9709
diff
changeset
|
61 return U; |
8869
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
62 } |
c3b743b1b1c6
preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents:
8367
diff
changeset
|
63 |
15039
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14854
diff
changeset
|
64 DEFUN (lu, args, nargout, |
21966
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
65 doc: /* -*- texinfo -*- |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
66 @deftypefn {} {[@var{L}, @var{U}] =} lu (@var{A}) |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
67 @deftypefnx {} {[@var{L}, @var{U}, @var{P}] =} lu (@var{A}) |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
68 @deftypefnx {} {[@var{L}, @var{U}, @var{P}, @var{Q}] =} lu (@var{S}) |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
69 @deftypefnx {} {[@var{L}, @var{U}, @var{P}, @var{Q}, @var{R}] =} lu (@var{S}) |
28154 | 70 @deftypefnx {} {[@dots{}] =} lu (@var{S}, @var{thresh}) |
21966
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
71 @deftypefnx {} {@var{y} =} lu (@dots{}) |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
72 @deftypefnx {} {[@dots{}] =} lu (@dots{}, "vector") |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
73 @cindex LU decomposition |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
74 Compute the LU@tie{}decomposition of @var{A}. |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
75 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
76 If @var{A} is full then subroutines from @sc{lapack} are used, and if |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
77 @var{A} is sparse then @sc{umfpack} is used. |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
78 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
79 The result is returned in a permuted form, according to the optional return |
28154 | 80 value @var{P}. For example, given the matrix @code{@var{A} = [1, 2; 3, 4]}, |
21966
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
81 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
82 @example |
28154 | 83 [@var{L}, @var{U}, @var{P}] = lu (@var{A}) |
21966
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
84 @end example |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
85 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
86 @noindent |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
87 returns |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
88 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
89 @example |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
90 @group |
28154 | 91 L = |
21966
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
92 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
93 1.00000 0.00000 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
94 0.33333 1.00000 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
95 |
28154 | 96 U = |
21966
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
97 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
98 3.00000 4.00000 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
99 0.00000 0.66667 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
100 |
28154 | 101 P = |
21966
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
102 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
103 0 1 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
104 1 0 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
105 @end group |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
106 @end example |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
107 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
108 The matrix is not required to be square. |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
109 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
110 When called with two or three output arguments and a sparse input matrix, |
25036
fa2f8ffd088e
Add note about the preferred way to call lu with sparse input (bug #53390).
Rik <rik@octave.org>
parents:
24534
diff
changeset
|
111 @code{lu} does not attempt to perform sparsity preserving column permutations. |
fa2f8ffd088e
Add note about the preferred way to call lu with sparse input (bug #53390).
Rik <rik@octave.org>
parents:
24534
diff
changeset
|
112 Called with a fourth output argument, the sparsity preserving column |
fa2f8ffd088e
Add note about the preferred way to call lu with sparse input (bug #53390).
Rik <rik@octave.org>
parents:
24534
diff
changeset
|
113 transformation @var{Q} is returned, such that |
fa2f8ffd088e
Add note about the preferred way to call lu with sparse input (bug #53390).
Rik <rik@octave.org>
parents:
24534
diff
changeset
|
114 @code{@var{P} * @var{A} * @var{Q} = @var{L} * @var{U}}. This is the |
fa2f8ffd088e
Add note about the preferred way to call lu with sparse input (bug #53390).
Rik <rik@octave.org>
parents:
24534
diff
changeset
|
115 @strong{preferred} way to call @code{lu} with sparse input matrices. |
21966
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
116 |
25036
fa2f8ffd088e
Add note about the preferred way to call lu with sparse input (bug #53390).
Rik <rik@octave.org>
parents:
24534
diff
changeset
|
117 Called with a fifth output argument and a sparse input matrix, @code{lu} |
fa2f8ffd088e
Add note about the preferred way to call lu with sparse input (bug #53390).
Rik <rik@octave.org>
parents:
24534
diff
changeset
|
118 attempts to use a scaling factor @var{R} on the input matrix such that |
21966
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
119 @code{@var{P} * (@var{R} \ @var{A}) * @var{Q} = @var{L} * @var{U}}. |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
120 This typically leads to a sparser and more stable factorization. |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
121 |
28154 | 122 An additional input argument @var{thresh} that defines the pivoting |
123 threshold can be given. @var{thresh} can be a scalar, in which case | |
21966
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
124 it defines the @sc{umfpack} pivoting tolerance for both symmetric and |
28154 | 125 unsymmetric cases. If @var{thresh} is a 2-element vector, then the first |
21966
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
126 element defines the pivoting tolerance for the unsymmetric @sc{umfpack} |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
127 pivoting strategy and the second for the symmetric strategy. By default, |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
128 the values defined by @code{spparms} are used ([0.1, 0.001]). |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
129 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
130 Given the string argument @qcode{"vector"}, @code{lu} returns the values |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
131 of @var{P} and @var{Q} as vector values, such that for full matrix, |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
132 @code{@var{A}(@var{P},:) = @var{L} * @var{U}}, and @code{@var{R}(@var{P},:) |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
133 * @var{A}(:,@var{Q}) = @var{L} * @var{U}}. |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
134 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
135 With two output arguments, returns the permuted forms of the upper and |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
136 lower triangular matrices, such that @code{@var{A} = @var{L} * @var{U}}. |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
137 With one output argument @var{y}, then the matrix returned by the |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
138 @sc{lapack} routines is returned. If the input matrix is sparse then the |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
139 matrix @var{L} is embedded into @var{U} to give a return value similar to |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
140 the full case. For both full and sparse matrices, @code{lu} loses the |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
141 permutation information. |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
142 @seealso{luupdate, ilu, chol, hess, qr, qz, schur, svd} |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
143 @end deftypefn */) |
2928 | 144 { |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
145 int nargin = args.length (); |
23583
b7747a2c88b2
maint: Deprecate is_sparse_type and replace with issparse.
Rik <rik@octave.org>
parents:
23582
diff
changeset
|
146 bool issparse = (nargin > 0 && args(0).issparse ()); |
2928 | 147 |
20909
03e4ddd49396
omit unnecessary nargout checks
John W. Eaton <jwe@octave.org>
parents:
20898
diff
changeset
|
148 if (nargin < 1 || (issparse && nargin > 3) || (! issparse && nargin > 2)) |
20802
8bb38ba1bad6
eliminate return statements after calls to print_usage
John W. Eaton <jwe@octave.org>
parents:
20711
diff
changeset
|
149 print_usage (); |
2928 | 150 |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
151 bool vecout = false; |
28154 | 152 Matrix thresh; |
20892 | 153 int n = 1; |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
154 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
155 while (n < nargin) |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
156 { |
18112
b560bac0fca2
maint: Don't use space between 'args' and '(' when doing indexing.
Rik <rik@octave.org>
parents:
18100
diff
changeset
|
157 if (args(n).is_string ()) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
158 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
159 std::string tmp = args(n++).string_value (); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
160 |
20892 | 161 if (tmp == "vector") |
19743
67f2c76f9f4d
Remove unnecessary checking of error_state after is_string validation.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
162 vecout = true; |
67f2c76f9f4d
Remove unnecessary checking of error_state after is_string validation.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
163 else |
67f2c76f9f4d
Remove unnecessary checking of error_state after is_string validation.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
164 error ("lu: unrecognized string argument"); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
165 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
166 else |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
167 { |
20892 | 168 if (! issparse) |
28154 | 169 error ("lu: can not define pivoting threshold THRESH for full matrices"); |
21055
5e00ed38a58b
maint: Replace if/error/else paradigm with just if/error.
Rik <rik@octave.org>
parents:
20941
diff
changeset
|
170 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
171 Matrix tmp = args(n++).matrix_value (); |
20892 | 172 if (tmp.numel () == 1) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
173 { |
30390
a61e1a0f6024
maint: style check C++ files in libinterp/ ahead of 7.1 release.
Rik <rik@octave.org>
parents:
29961
diff
changeset
|
174 thresh.resize (1, 2); |
28154 | 175 thresh(0) = tmp(0); |
176 thresh(1) = tmp(0); | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
177 } |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
178 else if (tmp.numel () == 2) |
28154 | 179 thresh = tmp; |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
180 else |
28154 | 181 error ("lu: THRESH must be a 1- or 2-element vector"); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
182 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
183 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
184 |
20892 | 185 octave_value_list retval; |
186 bool scale = (nargout == 5); | |
187 | |
2928 | 188 octave_value arg = args(0); |
189 | |
5275 | 190 octave_idx_type nr = arg.rows (); |
191 octave_idx_type nc = arg.columns (); | |
2928 | 192 |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
193 if (issparse) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
194 { |
23577
80c42f4cca13
maint: Deprecate is_empty and replace with isempty.
Rik <rik@octave.org>
parents:
23450
diff
changeset
|
195 if (arg.isempty ()) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
196 return octave_value_list (5, SparseMatrix ()); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
197 |
23582
0cc2011d800e
maint: Deprecate is_real_type and replace with isreal.
Rik <rik@octave.org>
parents:
23581
diff
changeset
|
198 if (arg.isreal ()) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
199 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
200 SparseMatrix m = arg.sparse_matrix_value (); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
201 |
18678
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
18497
diff
changeset
|
202 if (nargout < 4) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
203 { |
25037
1c077d652c57
Add new warning ID and message when lu is called with sparse input incorrectly (bug #53390).
Rik <rik@octave.org>
parents:
25036
diff
changeset
|
204 warning_with_id ("Octave:lu:sparse_input", |
1c077d652c57
Add new warning ID and message when lu is called with sparse input incorrectly (bug #53390).
Rik <rik@octave.org>
parents:
25036
diff
changeset
|
205 "lu: function may fail when called with less than 4 output arguments and a sparse input"); |
1c077d652c57
Add new warning ID and message when lu is called with sparse input incorrectly (bug #53390).
Rik <rik@octave.org>
parents:
25036
diff
changeset
|
206 |
20892 | 207 ColumnVector Qinit (nc); |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
208 for (octave_idx_type i = 0; i < nc; i++) |
20892 | 209 Qinit(i) = i; |
29961
7d6709900da7
eliminate octave:: namespace tags in DEFUN and DEFMETHOD and more
John W. Eaton <jwe@octave.org>
parents:
29958
diff
changeset
|
210 math::sparse_lu<SparseMatrix> fact (m, Qinit, thresh, false, |
7d6709900da7
eliminate octave:: namespace tags in DEFUN and DEFMETHOD and more
John W. Eaton <jwe@octave.org>
parents:
29958
diff
changeset
|
211 true); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
212 |
18678
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
18497
diff
changeset
|
213 if (nargout < 2) |
19861
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
214 retval(0) = fact.Y (); |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
215 else |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
216 { |
20921
4d3daf7e43f3
eliminate trailing whitespace in source files
John W. Eaton <jwe@octave.org>
parents:
20909
diff
changeset
|
217 retval.resize (nargout == 3 ? 3 : 2); |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
218 retval(1) |
27276
7455523fdf01
style fixes: avoid breaking lines immediately after '('
John W. Eaton <jwe@octave.org>
parents:
26376
diff
changeset
|
219 = octave_value (fact.U () * fact.Pc_mat ().transpose (), |
7455523fdf01
style fixes: avoid breaking lines immediately after '('
John W. Eaton <jwe@octave.org>
parents:
26376
diff
changeset
|
220 MatrixType (MatrixType::Permuted_Upper, |
7455523fdf01
style fixes: avoid breaking lines immediately after '('
John W. Eaton <jwe@octave.org>
parents:
26376
diff
changeset
|
221 nc, fact.col_perm ())); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
222 |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
223 PermMatrix P = fact.Pr_mat (); |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
224 SparseMatrix L = fact.L (); |
20892 | 225 |
226 if (nargout == 2) | |
22407
34ce5be04942
maint: Style check C++ code in libinterp/.
Rik <rik@octave.org>
parents:
22372
diff
changeset
|
227 retval(0) |
34ce5be04942
maint: Style check C++ code in libinterp/.
Rik <rik@octave.org>
parents:
22372
diff
changeset
|
228 = octave_value (P.transpose () * L, |
34ce5be04942
maint: Style check C++ code in libinterp/.
Rik <rik@octave.org>
parents:
22372
diff
changeset
|
229 MatrixType (MatrixType::Permuted_Lower, |
34ce5be04942
maint: Style check C++ code in libinterp/.
Rik <rik@octave.org>
parents:
22372
diff
changeset
|
230 nr, fact.row_perm ())); |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
231 else |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
232 { |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
233 retval(0) = L; |
18678
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
18497
diff
changeset
|
234 if (vecout) |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
235 retval(2) = fact.Pr_vec(); |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
236 else |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
237 retval(2) = P; |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
238 } |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
239 } |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
240 } |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
241 else |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
242 { |
20892 | 243 retval.resize (scale ? 5 : 4); |
29961
7d6709900da7
eliminate octave:: namespace tags in DEFUN and DEFMETHOD and more
John W. Eaton <jwe@octave.org>
parents:
29958
diff
changeset
|
244 math::sparse_lu<SparseMatrix> fact (m, thresh, scale); |
2928 | 245 |
20892 | 246 retval(0) = octave_value (fact.L (), |
247 MatrixType (MatrixType::Lower)); | |
248 retval(1) = octave_value (fact.U (), | |
249 MatrixType (MatrixType::Upper)); | |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
250 |
19861
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
251 if (vecout) |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
252 { |
20892 | 253 retval(2) = fact.Pr_vec (); |
19861
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
254 retval(3) = fact.Pc_vec (); |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
255 } |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
256 else |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
257 { |
20892 | 258 retval(2) = fact.Pr_mat (); |
19861
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
259 retval(3) = fact.Pc_mat (); |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
260 } |
20892 | 261 |
262 if (scale) | |
263 retval(4) = fact.R (); | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
264 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
265 } |
23581
c3075ae020e1
maint: Deprecate is_complex_type and replace with iscomplex.
Rik <rik@octave.org>
parents:
23577
diff
changeset
|
266 else if (arg.iscomplex ()) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
267 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
268 SparseComplexMatrix m = arg.sparse_complex_matrix_value (); |
2928 | 269 |
18678
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
18497
diff
changeset
|
270 if (nargout < 4) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
271 { |
25037
1c077d652c57
Add new warning ID and message when lu is called with sparse input incorrectly (bug #53390).
Rik <rik@octave.org>
parents:
25036
diff
changeset
|
272 warning_with_id ("Octave:lu:sparse_input", |
1c077d652c57
Add new warning ID and message when lu is called with sparse input incorrectly (bug #53390).
Rik <rik@octave.org>
parents:
25036
diff
changeset
|
273 "lu: function may fail when called with less than 4 output arguments and a sparse input"); |
1c077d652c57
Add new warning ID and message when lu is called with sparse input incorrectly (bug #53390).
Rik <rik@octave.org>
parents:
25036
diff
changeset
|
274 |
20892 | 275 ColumnVector Qinit (nc); |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
276 for (octave_idx_type i = 0; i < nc; i++) |
20892 | 277 Qinit(i) = i; |
29961
7d6709900da7
eliminate octave:: namespace tags in DEFUN and DEFMETHOD and more
John W. Eaton <jwe@octave.org>
parents:
29958
diff
changeset
|
278 math::sparse_lu<SparseComplexMatrix> fact (m, Qinit, thresh, |
31607
aac27ad79be6
maint: Re-indent code after switch to using namespace macros.
Rik <rik@octave.org>
parents:
31605
diff
changeset
|
279 false, true); |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
280 |
18678
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
18497
diff
changeset
|
281 if (nargout < 2) |
19861
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
282 retval(0) = fact.Y (); |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
283 else |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
284 { |
20921
4d3daf7e43f3
eliminate trailing whitespace in source files
John W. Eaton <jwe@octave.org>
parents:
20909
diff
changeset
|
285 retval.resize (nargout == 3 ? 3 : 2); |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
286 retval(1) |
27276
7455523fdf01
style fixes: avoid breaking lines immediately after '('
John W. Eaton <jwe@octave.org>
parents:
26376
diff
changeset
|
287 = octave_value (fact.U () * fact.Pc_mat ().transpose (), |
7455523fdf01
style fixes: avoid breaking lines immediately after '('
John W. Eaton <jwe@octave.org>
parents:
26376
diff
changeset
|
288 MatrixType (MatrixType::Permuted_Upper, |
7455523fdf01
style fixes: avoid breaking lines immediately after '('
John W. Eaton <jwe@octave.org>
parents:
26376
diff
changeset
|
289 nc, fact.col_perm ())); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
290 |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
291 PermMatrix P = fact.Pr_mat (); |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
292 SparseComplexMatrix L = fact.L (); |
20892 | 293 if (nargout == 2) |
19861
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
294 retval(0) |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
295 = octave_value (P.transpose () * L, |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
296 MatrixType (MatrixType::Permuted_Lower, |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
297 nr, fact.row_perm ())); |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
298 else |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
299 { |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
300 retval(0) = L; |
18678
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
18497
diff
changeset
|
301 if (vecout) |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
302 retval(2) = fact.Pr_vec(); |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
303 else |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
304 retval(2) = P; |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
305 } |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
306 } |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
307 } |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
308 else |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
309 { |
20892 | 310 retval.resize (scale ? 5 : 4); |
29961
7d6709900da7
eliminate octave:: namespace tags in DEFUN and DEFMETHOD and more
John W. Eaton <jwe@octave.org>
parents:
29958
diff
changeset
|
311 math::sparse_lu<SparseComplexMatrix> fact (m, thresh, scale); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
312 |
20892 | 313 retval(0) = octave_value (fact.L (), |
314 MatrixType (MatrixType::Lower)); | |
315 retval(1) = octave_value (fact.U (), | |
316 MatrixType (MatrixType::Upper)); | |
317 | |
318 if (vecout) | |
319 { | |
320 retval(2) = fact.Pr_vec (); | |
321 retval(3) = fact.Pc_vec (); | |
322 } | |
323 else | |
324 { | |
325 retval(2) = fact.Pr_mat (); | |
326 retval(3) = fact.Pc_mat (); | |
327 } | |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
328 |
19861
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
329 if (scale) |
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19743
diff
changeset
|
330 retval(4) = fact.R (); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
331 } |
18427
93c019f00e59
maint: whitespace fixes for libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
18426
diff
changeset
|
332 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
333 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
334 else |
21100
e39e05d90788
Switch gripe_XXX to either err_XXX or warn_XXX naming scheme.
Rik <rik@octave.org>
parents:
21055
diff
changeset
|
335 err_wrong_type_arg ("lu", arg); |
2928 | 336 } |
337 else | |
338 { | |
23577
80c42f4cca13
maint: Deprecate is_empty and replace with isempty.
Rik <rik@octave.org>
parents:
23450
diff
changeset
|
339 if (arg.isempty ()) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
340 return octave_value_list (3, Matrix ()); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
341 |
23582
0cc2011d800e
maint: Deprecate is_real_type and replace with isreal.
Rik <rik@octave.org>
parents:
23581
diff
changeset
|
342 if (arg.isreal ()) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
343 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
344 if (arg.is_single_type ()) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
345 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
346 FloatMatrix m = arg.float_matrix_value (); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
347 |
29961
7d6709900da7
eliminate octave:: namespace tags in DEFUN and DEFMETHOD and more
John W. Eaton <jwe@octave.org>
parents:
29958
diff
changeset
|
348 math::lu<FloatMatrix> fact (m); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
349 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
350 switch (nargout) |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
351 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
352 case 0: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
353 case 1: |
20892 | 354 retval = ovl (fact.Y ()); |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
355 break; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
356 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
357 case 2: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
358 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
359 PermMatrix P = fact.P (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
360 FloatMatrix L = P.transpose () * fact.L (); |
20892 | 361 retval = ovl (L, get_lu_u (fact)); |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
362 } |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
363 break; |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
364 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
365 case 3: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
366 default: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
367 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
368 if (vecout) |
20892 | 369 retval = ovl (get_lu_l (fact), get_lu_u (fact), |
370 fact.P_vec ()); | |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
371 else |
20892 | 372 retval = ovl (get_lu_l (fact), get_lu_u (fact), |
373 fact.P ()); | |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
374 } |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
375 break; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
376 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
377 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
378 else |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
379 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
380 Matrix m = arg.matrix_value (); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
381 |
29961
7d6709900da7
eliminate octave:: namespace tags in DEFUN and DEFMETHOD and more
John W. Eaton <jwe@octave.org>
parents:
29958
diff
changeset
|
382 math::lu<Matrix> fact (m); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
383 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
384 switch (nargout) |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
385 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
386 case 0: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
387 case 1: |
20892 | 388 retval = ovl (fact.Y ()); |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
389 break; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
390 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
391 case 2: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
392 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
393 PermMatrix P = fact.P (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
394 Matrix L = P.transpose () * fact.L (); |
20892 | 395 retval = ovl (L, get_lu_u (fact)); |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
396 } |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
397 break; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
398 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
399 case 3: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
400 default: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
401 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
402 if (vecout) |
20892 | 403 retval = ovl (get_lu_l (fact), get_lu_u (fact), |
404 fact.P_vec ()); | |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
405 else |
20892 | 406 retval = ovl (get_lu_l (fact), get_lu_u (fact), |
407 fact.P ()); | |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
408 } |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
409 break; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
410 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
411 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
412 } |
23581
c3075ae020e1
maint: Deprecate is_complex_type and replace with iscomplex.
Rik <rik@octave.org>
parents:
23577
diff
changeset
|
413 else if (arg.iscomplex ()) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
414 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
415 if (arg.is_single_type ()) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
416 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
417 FloatComplexMatrix m = arg.float_complex_matrix_value (); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
418 |
29961
7d6709900da7
eliminate octave:: namespace tags in DEFUN and DEFMETHOD and more
John W. Eaton <jwe@octave.org>
parents:
29958
diff
changeset
|
419 math::lu<FloatComplexMatrix> fact (m); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
420 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
421 switch (nargout) |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
422 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
423 case 0: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
424 case 1: |
20892 | 425 retval = ovl (fact.Y ()); |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
426 break; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
427 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
428 case 2: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
429 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
430 PermMatrix P = fact.P (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
431 FloatComplexMatrix L = P.transpose () * fact.L (); |
20892 | 432 retval = ovl (L, get_lu_u (fact)); |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
433 } |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
434 break; |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
435 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
436 case 3: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
437 default: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
438 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
439 if (vecout) |
20892 | 440 retval = ovl (get_lu_l (fact), get_lu_u (fact), |
441 fact.P_vec ()); | |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
442 else |
20892 | 443 retval = ovl (get_lu_l (fact), get_lu_u (fact), |
444 fact.P ()); | |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
445 } |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
446 break; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
447 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
448 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
449 else |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
450 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
451 ComplexMatrix m = arg.complex_matrix_value (); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
452 |
29961
7d6709900da7
eliminate octave:: namespace tags in DEFUN and DEFMETHOD and more
John W. Eaton <jwe@octave.org>
parents:
29958
diff
changeset
|
453 math::lu<ComplexMatrix> fact (m); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
454 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
455 switch (nargout) |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
456 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
457 case 0: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
458 case 1: |
20892 | 459 retval = ovl (fact.Y ()); |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
460 break; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
461 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
462 case 2: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
463 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
464 PermMatrix P = fact.P (); |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
465 ComplexMatrix L = P.transpose () * fact.L (); |
20892 | 466 retval = ovl (L, get_lu_u (fact)); |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
467 } |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
468 break; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7515
diff
changeset
|
469 |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
470 case 3: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
471 default: |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
472 { |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
473 if (vecout) |
20892 | 474 retval = ovl (get_lu_l (fact), get_lu_u (fact), |
475 fact.P_vec ()); | |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
476 else |
20892 | 477 retval = ovl (get_lu_l (fact), get_lu_u (fact), |
478 fact.P ()); | |
20555
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
479 } |
f90c8372b7ba
eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents:
20228
diff
changeset
|
480 break; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
481 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
482 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
10094
diff
changeset
|
483 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
484 else |
21100
e39e05d90788
Switch gripe_XXX to either err_XXX or warn_XXX naming scheme.
Rik <rik@octave.org>
parents:
21055
diff
changeset
|
485 err_wrong_type_arg ("lu", arg); |
2928 | 486 } |
487 | |
488 return retval; | |
489 } | |
490 | |
491 /* | |
21317
a4faec57f4c8
maint: remove semicolon after %!assert tests to follow Octave conventions.
Rik <rik@octave.org>
parents:
21301
diff
changeset
|
492 %!assert(lu ([1, 2; 3, 4]), [3, 4; 1/3, 2/3], eps) |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
493 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
494 %!test |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
495 %! [l, u] = lu ([1, 2; 3, 4]); |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
496 %! assert (l, [1/3, 1; 1, 0], sqrt (eps)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
497 %! assert (u, [3, 4; 0, 2/3], sqrt (eps)); |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
498 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
499 %!test |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
500 %! [l, u, p] = lu ([1, 2; 3, 4]); |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
501 %! assert (l, [1, 0; 1/3, 1], sqrt (eps)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
502 %! assert (u, [3, 4; 0, 2/3], sqrt (eps)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
503 %! assert (p(:,:), [0, 1; 1, 0], sqrt (eps)); |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
504 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
505 %!test |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
506 %! [l, u, p] = lu ([1, 2; 3, 4], "vector"); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
507 %! assert (l, [1, 0; 1/3, 1], sqrt (eps)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
508 %! assert (u, [3, 4; 0, 2/3], sqrt (eps)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
509 %! assert (p, [2;1], sqrt (eps)); |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
510 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
511 %!test |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
512 %! [l, u, p] = lu ([1, 2; 3, 4; 5, 6]); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
513 %! assert (l, [1, 0; 1/5, 1; 3/5, 1/2], sqrt (eps)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
514 %! assert (u, [5, 6; 0, 4/5], sqrt (eps)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
515 %! assert (p(:,:), [0, 0, 1; 1, 0, 0; 0 1 0], sqrt (eps)); |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
516 |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
517 %!assert (lu (single ([1, 2; 3, 4])), single ([3, 4; 1/3, 2/3]), eps ("single")) |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
518 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
519 %!test |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
520 %! [l, u] = lu (single ([1, 2; 3, 4])); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
521 %! assert (l, single ([1/3, 1; 1, 0]), sqrt (eps ("single"))); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
522 %! assert (u, single ([3, 4; 0, 2/3]), sqrt (eps ("single"))); |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
523 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
524 %!test |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
525 %! [l, u, p] = lu (single ([1, 2; 3, 4])); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
526 %! assert (l, single ([1, 0; 1/3, 1]), sqrt (eps ("single"))); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
527 %! assert (u, single ([3, 4; 0, 2/3]), sqrt (eps ("single"))); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
528 %! assert (p(:,:), single ([0, 1; 1, 0]), sqrt (eps ("single"))); |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
529 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
530 %!test |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
531 %! [l, u, p] = lu (single ([1, 2; 3, 4]), "vector"); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
532 %! assert (l, single ([1, 0; 1/3, 1]), sqrt (eps ("single"))); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
533 %! assert (u, single ([3, 4; 0, 2/3]), sqrt (eps ("single"))); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
534 %! assert (p, single ([2;1]), sqrt (eps ("single"))); |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
535 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
536 %!test |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
537 %! [l u p] = lu (single ([1, 2; 3, 4; 5, 6])); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
538 %! assert (l, single ([1, 0; 1/5, 1; 3/5, 1/2]), sqrt (eps ("single"))); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
539 %! assert (u, single ([5, 6; 0, 4/5]), sqrt (eps ("single"))); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
540 %! assert (p(:,:), single ([0, 0, 1; 1, 0, 0; 0 1 0]), sqrt (eps ("single"))); |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
541 |
18488
4c2465444a96
Use %!testif HAVE_UMFPACK for sparse lu test added in cset 2a45b6b87bee
Mike Miller <mtmiller@ieee.org>
parents:
18427
diff
changeset
|
542 %!testif HAVE_UMFPACK |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
543 %! Bi = [1 2 3 4 5 2 3 6 7 8 4 5 7 8 9]; |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
544 %! Bj = [1 3 4 5 6 7 8 9 11 12 13 14 15 16 17]; |
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
545 %! Bv = [1 1 1 1 1 1 -1 1 1 1 1 -1 1 -1 1]; |
18427
93c019f00e59
maint: whitespace fixes for libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
18426
diff
changeset
|
546 %! B = sparse (Bi, Bj, Bv); |
25037
1c077d652c57
Add new warning ID and message when lu is called with sparse input incorrectly (bug #53390).
Rik <rik@octave.org>
parents:
25036
diff
changeset
|
547 %! warning ("off", "Octave:lu:sparse_input", "local"); |
18427
93c019f00e59
maint: whitespace fixes for libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
18426
diff
changeset
|
548 %! [L, U] = lu (B); |
18426
807e21ebddd5
maint: Style fixes for tests in libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
18425
diff
changeset
|
549 %! assert (L*U, B); |
18427
93c019f00e59
maint: whitespace fixes for libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
18426
diff
changeset
|
550 %! [L, U, P] = lu(B); |
18426
807e21ebddd5
maint: Style fixes for tests in libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
18425
diff
changeset
|
551 %! assert (P'*L*U, B); |
18427
93c019f00e59
maint: whitespace fixes for libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
18426
diff
changeset
|
552 %! [L, U, P, Q] = lu (B); |
18426
807e21ebddd5
maint: Style fixes for tests in libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
18425
diff
changeset
|
553 %! assert (P'*L*U*Q', B); |
18425
2a45b6b87bee
correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents:
17787
diff
changeset
|
554 |
25037
1c077d652c57
Add new warning ID and message when lu is called with sparse input incorrectly (bug #53390).
Rik <rik@octave.org>
parents:
25036
diff
changeset
|
555 %!error lu () |
25095
71ed409b2ffb
Change lu BIST test to use testif HAVE_UMFPACK.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
556 %!testif HAVE_UMFPACK |
71ed409b2ffb
Change lu BIST test to use testif HAVE_UMFPACK.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
557 %! fail ("[l,u] = lu (sparse (magic (3)))", "warning", "function may fail"); |
25037
1c077d652c57
Add new warning ID and message when lu is called with sparse input incorrectly (bug #53390).
Rik <rik@octave.org>
parents:
25036
diff
changeset
|
558 %!error <can not define pivoting threshold> lu ([1, 2; 3, 4], 2) |
1c077d652c57
Add new warning ID and message when lu is called with sparse input incorrectly (bug #53390).
Rik <rik@octave.org>
parents:
25036
diff
changeset
|
559 |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
560 */ |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
561 |
9708 | 562 static |
563 bool check_lu_dims (const octave_value& l, const octave_value& u, | |
564 const octave_value& p) | |
565 { | |
18100
6a71e5030df5
Follow coding convention of defining and initializing only 1 variable per line in liboctinterp.
Rik <rik@octave.org>
parents:
17787
diff
changeset
|
566 octave_idx_type m = l.rows (); |
6a71e5030df5
Follow coding convention of defining and initializing only 1 variable per line in liboctinterp.
Rik <rik@octave.org>
parents:
17787
diff
changeset
|
567 octave_idx_type k = u.rows (); |
6a71e5030df5
Follow coding convention of defining and initializing only 1 variable per line in liboctinterp.
Rik <rik@octave.org>
parents:
17787
diff
changeset
|
568 octave_idx_type n = u.columns (); |
20892 | 569 |
9708 | 570 return ((l.ndims () == 2 && u.ndims () == 2 && k == l.columns ()) |
19864
17d647821d61
maint: More cleanup of C++ code to follow Octave coding conventions.
John W. Eaton <jwe@octave.org>
parents:
19861
diff
changeset
|
571 && k == std::min (m, n) |
17d647821d61
maint: More cleanup of C++ code to follow Octave coding conventions.
John W. Eaton <jwe@octave.org>
parents:
19861
diff
changeset
|
572 && (p.is_undefined () || p.rows () == m)); |
9708 | 573 } |
574 | |
15039
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14854
diff
changeset
|
575 DEFUN (luupdate, args, , |
21966
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
576 doc: /* -*- texinfo -*- |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
577 @deftypefn {} {[@var{L}, @var{U}] =} luupdate (@var{L}, @var{U}, @var{x}, @var{y}) |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
578 @deftypefnx {} {[@var{L}, @var{U}, @var{P}] =} luupdate (@var{L}, @var{U}, @var{P}, @var{x}, @var{y}) |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
579 Given an LU@tie{}factorization of a real or complex matrix |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
580 @w{@var{A} = @var{L}*@var{U}}, @var{L}@tie{}lower unit trapezoidal and |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
581 @var{U}@tie{}upper trapezoidal, return the LU@tie{}factorization |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
582 of @w{@var{A} + @var{x}*@var{y}.'}, where @var{x} and @var{y} are |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
583 column vectors (rank-1 update) or matrices with equal number of columns |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
584 (rank-k update). |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
585 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
586 Optionally, row-pivoted updating can be used by supplying a row permutation |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
587 (pivoting) matrix @var{P}; in that case, an updated permutation matrix is |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
588 returned. Note that if @var{L}, @var{U}, @var{P} is a pivoted |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
589 LU@tie{}factorization as obtained by @code{lu}: |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
590 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
591 @example |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
592 [@var{L}, @var{U}, @var{P}] = lu (@var{A}); |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
593 @end example |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
594 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
595 @noindent |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
596 then a factorization of @tcode{@var{A}+@var{x}*@var{y}.'} can be obtained |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
597 either as |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
598 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
599 @example |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
600 [@var{L1}, @var{U1}] = lu (@var{L}, @var{U}, @var{P}*@var{x}, @var{y}) |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
601 @end example |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
602 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
603 @noindent |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
604 or |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
605 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
606 @example |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
607 [@var{L1}, @var{U1}, @var{P1}] = lu (@var{L}, @var{U}, @var{P}, @var{x}, @var{y}) |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
608 @end example |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
609 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
610 The first form uses the unpivoted algorithm, which is faster, but less |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
611 stable. The second form uses a slower pivoted algorithm, which is more |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
612 stable. |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
613 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
614 The matrix case is done as a sequence of rank-1 updates; thus, for large |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
615 enough k, it will be both faster and more accurate to recompute the |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
616 factorization from scratch. |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
617 @seealso{lu, cholupdate, qrupdate} |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21724
diff
changeset
|
618 @end deftypefn */) |
9708 | 619 { |
20816
b16bcd7a2a33
Use int rather than octave_idx_type for nargin data type.
Rik <rik@octave.org>
parents:
20812
diff
changeset
|
620 int nargin = args.length (); |
9708 | 621 |
622 if (nargin != 4 && nargin != 5) | |
20802
8bb38ba1bad6
eliminate return statements after calls to print_usage
John W. Eaton <jwe@octave.org>
parents:
20711
diff
changeset
|
623 print_usage (); |
9708 | 624 |
20812
d9ca869ca124
maint: Clean-up more instances of print_usage().
Rik <rik@octave.org>
parents:
20802
diff
changeset
|
625 bool pivoted = (nargin == 5); |
d9ca869ca124
maint: Clean-up more instances of print_usage().
Rik <rik@octave.org>
parents:
20802
diff
changeset
|
626 |
9708 | 627 octave_value argl = args(0); |
628 octave_value argu = args(1); | |
23450
855122b993da
maint: Wrap tertiary operator in parentheses "(COND ? x : y)".
Rik <rik@octave.org>
parents:
23220
diff
changeset
|
629 octave_value argp = (pivoted ? args(2) : octave_value ()); |
9708 | 630 octave_value argx = args(2 + pivoted); |
631 octave_value argy = args(3 + pivoted); | |
632 | |
23586
f6c5db0a02e7
maint: Deprecate is_numeric_type and replace with isnumeric.
Rik <rik@octave.org>
parents:
23583
diff
changeset
|
633 if (! (argl.isnumeric () && argu.isnumeric () |
f6c5db0a02e7
maint: Deprecate is_numeric_type and replace with isnumeric.
Rik <rik@octave.org>
parents:
23583
diff
changeset
|
634 && argx.isnumeric () && argy.isnumeric () |
20892 | 635 && (! pivoted || argp.is_perm_matrix ()))) |
636 error ("luupdate: L, U, X, and Y must be numeric"); | |
9708 | 637 |
20892 | 638 if (! check_lu_dims (argl, argu, argp)) |
639 error ("luupdate: dimension mismatch"); | |
9708 | 640 |
20892 | 641 PermMatrix P = (pivoted |
642 ? argp.perm_matrix_value () | |
643 : PermMatrix::eye (argl.rows ())); | |
9708 | 644 |
23582
0cc2011d800e
maint: Deprecate is_real_type and replace with isreal.
Rik <rik@octave.org>
parents:
23581
diff
changeset
|
645 if (argl.isreal () && argu.isreal () |
0cc2011d800e
maint: Deprecate is_real_type and replace with isreal.
Rik <rik@octave.org>
parents:
23581
diff
changeset
|
646 && argx.isreal () && argy.isreal ()) |
20892 | 647 { |
648 // all real case | |
649 if (argl.is_single_type () || argu.is_single_type () | |
650 || argx.is_single_type () || argy.is_single_type ()) | |
651 { | |
652 FloatMatrix L = argl.float_matrix_value (); | |
653 FloatMatrix U = argu.float_matrix_value (); | |
654 FloatMatrix x = argx.float_matrix_value (); | |
655 FloatMatrix y = argy.float_matrix_value (); | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
656 |
29961
7d6709900da7
eliminate octave:: namespace tags in DEFUN and DEFMETHOD and more
John W. Eaton <jwe@octave.org>
parents:
29958
diff
changeset
|
657 math::lu<FloatMatrix> fact (L, U, P); |
20892 | 658 if (pivoted) |
659 fact.update_piv (x, y); | |
660 else | |
661 fact.update (x, y); | |
9708 | 662 |
20892 | 663 if (pivoted) |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20892
diff
changeset
|
664 return ovl (get_lu_l (fact), get_lu_u (fact), fact.P ()); |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20892
diff
changeset
|
665 else |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20892
diff
changeset
|
666 return ovl (get_lu_l (fact), get_lu_u (fact)); |
9708 | 667 } |
668 else | |
20892 | 669 { |
670 Matrix L = argl.matrix_value (); | |
671 Matrix U = argu.matrix_value (); | |
672 Matrix x = argx.matrix_value (); | |
673 Matrix y = argy.matrix_value (); | |
674 | |
29961
7d6709900da7
eliminate octave:: namespace tags in DEFUN and DEFMETHOD and more
John W. Eaton <jwe@octave.org>
parents:
29958
diff
changeset
|
675 math::lu<Matrix> fact (L, U, P); |
20892 | 676 if (pivoted) |
677 fact.update_piv (x, y); | |
678 else | |
679 fact.update (x, y); | |
680 | |
681 if (pivoted) | |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20892
diff
changeset
|
682 return ovl (get_lu_l (fact), get_lu_u (fact), fact.P ()); |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20892
diff
changeset
|
683 else |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20892
diff
changeset
|
684 return ovl (get_lu_l (fact), get_lu_u (fact)); |
20892 | 685 } |
9708 | 686 } |
687 else | |
20892 | 688 { |
689 // complex case | |
690 if (argl.is_single_type () || argu.is_single_type () | |
691 || argx.is_single_type () || argy.is_single_type ()) | |
692 { | |
693 FloatComplexMatrix L = argl.float_complex_matrix_value (); | |
694 FloatComplexMatrix U = argu.float_complex_matrix_value (); | |
695 FloatComplexMatrix x = argx.float_complex_matrix_value (); | |
696 FloatComplexMatrix y = argy.float_complex_matrix_value (); | |
697 | |
29961
7d6709900da7
eliminate octave:: namespace tags in DEFUN and DEFMETHOD and more
John W. Eaton <jwe@octave.org>
parents:
29958
diff
changeset
|
698 math::lu<FloatComplexMatrix> fact (L, U, P); |
20892 | 699 if (pivoted) |
700 fact.update_piv (x, y); | |
701 else | |
702 fact.update (x, y); | |
703 | |
704 if (pivoted) | |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20892
diff
changeset
|
705 return ovl (get_lu_l (fact), get_lu_u (fact), fact.P ()); |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20892
diff
changeset
|
706 else |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20892
diff
changeset
|
707 return ovl (get_lu_l (fact), get_lu_u (fact)); |
20892 | 708 } |
709 else | |
710 { | |
711 ComplexMatrix L = argl.complex_matrix_value (); | |
712 ComplexMatrix U = argu.complex_matrix_value (); | |
713 ComplexMatrix x = argx.complex_matrix_value (); | |
714 ComplexMatrix y = argy.complex_matrix_value (); | |
715 | |
29961
7d6709900da7
eliminate octave:: namespace tags in DEFUN and DEFMETHOD and more
John W. Eaton <jwe@octave.org>
parents:
29958
diff
changeset
|
716 math::lu<ComplexMatrix> fact (L, U, P); |
20892 | 717 if (pivoted) |
718 fact.update_piv (x, y); | |
719 else | |
720 fact.update (x, y); | |
721 | |
722 if (pivoted) | |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20892
diff
changeset
|
723 return ovl (get_lu_l (fact), get_lu_u (fact), fact.P ()); |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20892
diff
changeset
|
724 else |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20892
diff
changeset
|
725 return ovl (get_lu_l (fact), get_lu_u (fact)); |
20892 | 726 } |
727 } | |
9708 | 728 } |
729 | |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
730 /* |
10080
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
731 %!shared A, u, v, Ac, uc, vc |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
732 %! A = [0.091364 0.613038 0.999083; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
733 %! 0.594638 0.425302 0.603537; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
734 %! 0.383594 0.291238 0.085574; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
735 %! 0.265712 0.268003 0.238409; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
736 %! 0.669966 0.743851 0.445057 ]; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
737 %! |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
738 %! u = [0.85082; |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
739 %! 0.76426; |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
740 %! 0.42883; |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
741 %! 0.53010; |
10080
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
742 %! 0.80683 ]; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
743 %! |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
744 %! v = [0.98810; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
745 %! 0.24295; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
746 %! 0.43167 ]; |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
747 %! |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
748 %! Ac = [0.620405 + 0.956953i 0.480013 + 0.048806i 0.402627 + 0.338171i; |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
749 %! 0.589077 + 0.658457i 0.013205 + 0.279323i 0.229284 + 0.721929i; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
750 %! 0.092758 + 0.345687i 0.928679 + 0.241052i 0.764536 + 0.832406i; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
751 %! 0.912098 + 0.721024i 0.049018 + 0.269452i 0.730029 + 0.796517i; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
752 %! 0.112849 + 0.603871i 0.486352 + 0.142337i 0.355646 + 0.151496i ]; |
10080
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
753 %! |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
754 %! uc = [0.20351 + 0.05401i; |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
755 %! 0.13141 + 0.43708i; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
756 %! 0.29808 + 0.08789i; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
757 %! 0.69821 + 0.38844i; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
758 %! 0.74871 + 0.25821i ]; |
10080
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
759 %! |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
760 %! vc = [0.85839 + 0.29468i; |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
761 %! 0.20820 + 0.93090i; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
762 %! 0.86184 + 0.34689i ]; |
10080
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
763 %! |
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
764 |
10094
ab1101011a6d
lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents:
10080
diff
changeset
|
765 %!testif HAVE_QRUPDATE_LUU |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
766 %! [L,U,P] = lu (A); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
767 %! [L,U] = luupdate (L,U,P*u,v); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
768 %! assert (norm (vec (tril (L)-L), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
769 %! assert (norm (vec (triu (U)-U), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
770 %! assert (norm (vec (P'*L*U - A - u*v.'), Inf) < norm (A)*1e1*eps); |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
771 %! |
10094
ab1101011a6d
lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents:
10080
diff
changeset
|
772 %!testif HAVE_QRUPDATE_LUU |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
773 %! [L,U,P] = lu (Ac); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
774 %! [L,U] = luupdate (L,U,P*uc,vc); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
775 %! assert (norm (vec (tril (L)-L), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
776 %! assert (norm (vec (triu (U)-U), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
777 %! assert (norm (vec (P'*L*U - Ac - uc*vc.'), Inf) < norm (Ac)*1e1*eps); |
10080
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
778 |
10094
ab1101011a6d
lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents:
10080
diff
changeset
|
779 %!testif HAVE_QRUPDATE_LUU |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
780 %! [L,U,P] = lu (single (A)); |
31607
aac27ad79be6
maint: Re-indent code after switch to using namespace macros.
Rik <rik@octave.org>
parents:
31605
diff
changeset
|
781 %! [L,U] = luupdate (L,U,P* single (u), single (v)); |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
782 %! assert (norm (vec (tril (L)-L), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
783 %! assert (norm (vec (triu (U)-U), Inf) == 0); |
31607
aac27ad79be6
maint: Re-indent code after switch to using namespace macros.
Rik <rik@octave.org>
parents:
31605
diff
changeset
|
784 %! assert (norm (vec (P'*L*U - single (A) - single (u)* single (v).'), Inf) |
aac27ad79be6
maint: Re-indent code after switch to using namespace macros.
Rik <rik@octave.org>
parents:
31605
diff
changeset
|
785 %! < norm (single (A))*1e1* eps ("single")); |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
786 %! |
10094
ab1101011a6d
lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents:
10080
diff
changeset
|
787 %!testif HAVE_QRUPDATE_LUU |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
788 %! [L,U,P] = lu (single (Ac)); |
31607
aac27ad79be6
maint: Re-indent code after switch to using namespace macros.
Rik <rik@octave.org>
parents:
31605
diff
changeset
|
789 %! [L,U] = luupdate (L,U,P* single (uc),single (vc)); |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
790 %! assert (norm (vec (tril (L)-L), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
791 %! assert (norm (vec (triu (U)-U), Inf) == 0); |
31607
aac27ad79be6
maint: Re-indent code after switch to using namespace macros.
Rik <rik@octave.org>
parents:
31605
diff
changeset
|
792 %! assert (norm (vec (P'*L*U - single (Ac) - single (uc)* single (vc).'), Inf) |
aac27ad79be6
maint: Re-indent code after switch to using namespace macros.
Rik <rik@octave.org>
parents:
31605
diff
changeset
|
793 %! < norm (single (Ac))*1e1* eps ("single")); |
10080
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
794 |
10094
ab1101011a6d
lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents:
10080
diff
changeset
|
795 %!testif HAVE_QRUPDATE_LUU |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
796 %! [L,U,P] = lu (A); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
797 %! [L,U,P] = luupdate (L,U,P,u,v); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
798 %! assert (norm (vec (tril (L)-L), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
799 %! assert (norm (vec (triu (U)-U), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
800 %! assert (norm (vec (P'*L*U - A - u*v.'), Inf) < norm (A)*1e1*eps); |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
801 %! |
10094
ab1101011a6d
lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents:
10080
diff
changeset
|
802 %!testif HAVE_QRUPDATE_LUU |
18849
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18678
diff
changeset
|
803 %! [L,U,P] = lu (A); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18678
diff
changeset
|
804 %! [~,ordcols] = max (P,[],1); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18678
diff
changeset
|
805 %! [~,ordrows] = max (P,[],2); |
28915
c40a367a84c0
maint: Use Octave convention of space after function name in libinterp/.
Rik <rik@octave.org>
parents:
28888
diff
changeset
|
806 %! P1 = eye (size (P))(:,ordcols); |
c40a367a84c0
maint: Use Octave convention of space after function name in libinterp/.
Rik <rik@octave.org>
parents:
28888
diff
changeset
|
807 %! P2 = eye (size (P))(ordrows,:); |
c40a367a84c0
maint: Use Octave convention of space after function name in libinterp/.
Rik <rik@octave.org>
parents:
28888
diff
changeset
|
808 %! assert (P1 == P); |
c40a367a84c0
maint: Use Octave convention of space after function name in libinterp/.
Rik <rik@octave.org>
parents:
28888
diff
changeset
|
809 %! assert (P2 == P); |
18849
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18678
diff
changeset
|
810 %! [L,U,P] = luupdate (L,U,P,u,v); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18678
diff
changeset
|
811 %! [L,U,P1] = luupdate (L,U,P1,u,v); |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18678
diff
changeset
|
812 %! [L,U,P2] = luupdate (L,U,P2,u,v); |
28915
c40a367a84c0
maint: Use Octave convention of space after function name in libinterp/.
Rik <rik@octave.org>
parents:
28888
diff
changeset
|
813 %! assert (P1 == P); |
c40a367a84c0
maint: Use Octave convention of space after function name in libinterp/.
Rik <rik@octave.org>
parents:
28888
diff
changeset
|
814 %! assert (P2 == P); |
18849
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18678
diff
changeset
|
815 %! |
aa9ca67f09fb
make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents:
18678
diff
changeset
|
816 %!testif HAVE_QRUPDATE_LUU |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
817 %! [L,U,P] = lu (Ac); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
818 %! [L,U,P] = luupdate (L,U,P,uc,vc); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
819 %! assert (norm (vec (tril (L)-L), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
820 %! assert (norm (vec (triu (U)-U), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
821 %! assert (norm (vec (P'*L*U - Ac - uc*vc.'), Inf) < norm (Ac)*1e1*eps); |
10080
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
822 |
10094
ab1101011a6d
lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents:
10080
diff
changeset
|
823 %!testif HAVE_QRUPDATE_LUU |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
824 %! [L,U,P] = lu (single (A)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
825 %! [L,U,P] = luupdate (L,U,P,single (u),single (v)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
826 %! assert (norm (vec (tril (L)-L), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
827 %! assert (norm (vec (triu (U)-U), Inf) == 0); |
31607
aac27ad79be6
maint: Re-indent code after switch to using namespace macros.
Rik <rik@octave.org>
parents:
31605
diff
changeset
|
828 %! assert (norm (vec (P'*L*U - single (A) - single (u)* single (v).'), Inf) |
aac27ad79be6
maint: Re-indent code after switch to using namespace macros.
Rik <rik@octave.org>
parents:
31605
diff
changeset
|
829 %! < norm (single (A))*1e1* eps ("single")); |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11572
diff
changeset
|
830 %! |
10094
ab1101011a6d
lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents:
10080
diff
changeset
|
831 %!testif HAVE_QRUPDATE_LUU |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
832 %! [L,U,P] = lu (single (Ac)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
833 %! [L,U,P] = luupdate (L,U,P,single (uc),single (vc)); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
834 %! assert (norm (vec (tril (L)-L), Inf) == 0); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14361
diff
changeset
|
835 %! assert (norm (vec (triu (U)-U), Inf) == 0); |
31607
aac27ad79be6
maint: Re-indent code after switch to using namespace macros.
Rik <rik@octave.org>
parents:
31605
diff
changeset
|
836 %! assert (norm (vec (P'*L*U - single (Ac) - single (uc)* single (vc).'), Inf) |
aac27ad79be6
maint: Re-indent code after switch to using namespace macros.
Rik <rik@octave.org>
parents:
31605
diff
changeset
|
837 %! < norm (single (Ac))*1e1* eps ("single")); |
10080
cf70ee43077c
add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents:
9724
diff
changeset
|
838 */ |
29958
32c3a5805893
move DEFUN and DEFMETHOD functions inside octave namespace
John W. Eaton <jwe@octave.org>
parents:
29359
diff
changeset
|
839 |
31605
e88a07dec498
maint: Use macros to begin/end C++ namespaces.
Rik <rik@octave.org>
parents:
30565
diff
changeset
|
840 OCTAVE_END_NAMESPACE(octave) |