annotate libinterp/corefcn/lu.cc @ 20940:48b2ad5ee801

maint: Rename oct-obj.[cc|h] to ovl.[cc|h] for clarity. * libinterp/corefcn/oct-obj.h: Replaced. Contains a #warning issued if used in compilation and includes "ovl.h" instead. Original file renamed to libinterp/octave-value/ovl.h * libinterp/corefcn/oct-obj.cc: Renamed to libinterp/octave-value/ovl.cc. * oct-obj.cc (ovl ()): Added new function to return empty octave_value_list. * libinterp/corefcn/module.mk: Remove oct-obj.cc from build system. * libinterp/octave-value/module.mk: Add ovl.h and ovl.cc to build system. * mk-opts.pl, annotation-dialog.h, Cell.cc, __contourc__.cc, __dsearchn__.cc, __lin_interpn__.cc, __pchip_deriv__.cc, __qp__.cc, balance.cc, besselj.cc, betainc.cc, colloc.cc, conv2.cc, daspk.cc, dasrt.cc, dassl.cc, data.cc, debug.cc, defaults.cc, defun.cc, det.cc, dirfns.cc, dlmread.cc, eig.cc, error.cc, fft.cc, fft2.cc, fftn.cc, file-io.cc, filter.cc, find.cc, gammainc.cc, gcd.cc, getgrent.cc, getpwent.cc, getrusage.cc, givens.cc, graphics.cc, gripes.cc, help.cc, hess.cc, hex2num.cc, hook-fcn.h, input.cc, input.h, inv.cc, kron.cc, load-save.cc, lookup.cc, ls-hdf5.cc, ls-mat-ascii.cc, ls-mat4.cc, ls-mat5.cc, ls-oct-binary.cc, ls-oct-text.cc, lsode.cc, lu.cc, luinc.cc, max.cc, mex.cc, oct-hist.cc, oct-lvalue.cc, oct-lvalue.h, oct-map.h, oct-stream.cc, oct.h, ordschur.cc, pager.cc, pinv.cc, pr-output.cc, quad.cc, quadcc.cc, qz.cc, rand.cc, rcond.cc, regexp.cc, schur.cc, sparse-xpow.cc, strfns.cc, sub2ind.cc, svd.cc, sylvester.cc, symtab.h, syscalls.cc, sysdep.cc, time.cc, toplev.cc, tril.cc, tsearch.cc, typecast.cc, urlwrite.cc, utils.cc, variables.cc, xpow.cc, __delaunayn__.cc, __glpk__.cc, __voronoi__.cc, audiodevinfo.cc, audioread.cc, chol.cc, convhulln.cc, dmperm.cc, qr.cc, symbfact.cc, mkbuiltins, ov-base-diag.h, ov-base-int.cc, ov-base-mat.cc, ov-base-mat.h, ov-base-scalar.cc, ov-base-sparse.cc, ov-base-sparse.h, ov-base.cc, ov-bool-mat.cc, ov-bool.cc, ov-builtin.cc, ov-cell.cc, ov-colon.cc, ov-complex.cc, ov-cs-list.h, ov-cx-mat.cc, ov-dld-fcn.cc, ov-fcn.cc, ov-fcn.h, ov-float.cc, ov-flt-complex.cc, ov-flt-cx-mat.cc, ov-flt-re-mat.cc, ov-int16.cc, ov-int32.cc, ov-int64.cc, ov-int8.cc, ov-java.h, ov-mex-fcn.cc, ov-perm.h, ov-range.cc, ov-re-mat.cc, ov-scalar.cc, ov-str-mat.cc, ov-uint16.cc, ov-uint32.cc, ov-uint64.cc, ov-uint8.cc, ov-usr-fcn.cc, ov-usr-fcn.h, ov.cc, octave.cc, op-b-b.cc, op-b-bm.cc, op-b-sbm.cc, op-bm-b.cc, op-bm-bm.cc, op-bm-sbm.cc, op-cdm-cdm.cc, op-cell.cc, op-chm.cc, op-class.cc, op-cm-cm.cc, op-cm-cs.cc, op-cm-m.cc, op-cm-s.cc, op-cm-scm.cc, op-cm-sm.cc, op-cs-cm.cc, op-cs-cs.cc, op-cs-m.cc, op-cs-s.cc, op-cs-scm.cc, op-cs-sm.cc, op-dm-dm.cc, op-dm-scm.cc, op-dm-sm.cc, op-double-conv.cc, op-fcdm-fcdm.cc, op-fcdm-fdm.cc, op-fcm-fcm.cc, op-fcm-fcs.cc, op-fcm-fm.cc, op-fcm-fs.cc, op-fcn.cc, op-fcs-fcm.cc, op-fcs-fcs.cc, op-fcs-fm.cc, op-fcs-fs.cc, op-fdm-fdm.cc, op-float-conv.cc, op-fm-fcm.cc, op-fm-fcs.cc, op-fm-fm.cc, op-fm-fs.cc, op-fs-fcm.cc, op-fs-fcs.cc, op-fs-fm.cc, op-fs-fs.cc, op-i16-i16.cc, op-i32-i32.cc, op-i64-i64.cc, op-i8-i8.cc, op-int-concat.cc, op-int-conv.cc, op-m-cm.cc, op-m-cs.cc, op-m-m.cc, op-m-s.cc, op-m-scm.cc, op-m-sm.cc, op-pm-pm.cc, op-pm-scm.cc, op-pm-sm.cc, op-range.cc, op-s-cm.cc, op-s-cs.cc, op-s-m.cc, op-s-s.cc, op-s-scm.cc, op-s-sm.cc, op-sbm-b.cc, op-sbm-bm.cc, op-sbm-sbm.cc, op-scm-cm.cc, op-scm-cs.cc, op-scm-m.cc, op-scm-s.cc, op-scm-scm.cc, op-scm-sm.cc, op-sm-cm.cc, op-sm-cs.cc, op-sm-m.cc, op-sm-s.cc, op-sm-scm.cc, op-sm-sm.cc, op-str-m.cc, op-str-s.cc, op-str-str.cc, op-struct.cc, op-ui16-ui16.cc, op-ui32-ui32.cc, op-ui64-ui64.cc, op-ui8-ui8.cc, parse.h, pt-arg-list.cc, pt-assign.cc, pt-binop.cc, pt-cbinop.cc, pt-cell.cc, pt-colon.cc, pt-const.cc, pt-eval.h, pt-fcn-handle.cc, pt-funcall.h, pt-id.cc, pt-idx.cc, pt-jump.cc, pt-mat.cc, pt-select.cc, pt-unop.cc, token.cc, Array-sym.cc: replace '#include "oct-obj.h"' with '#include "ovl.h"'.
author Rik <rik@octave.org>
date Fri, 18 Dec 2015 16:04:56 -0800
parents 4d3daf7e43f3
children a4f5da7c5463
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
2928
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
1 /*
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
2
19697
4197fc428c7d maint: Update copyright notices for 2015.
John W. Eaton <jwe@octave.org>
parents: 19596
diff changeset
3 Copyright (C) 1996-2015 John W. Eaton
2928
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
4
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
5 This file is part of Octave.
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
6
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
7 Octave is free software; you can redistribute it and/or modify it
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
8 under the terms of the GNU General Public License as published by the
7016
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 5823
diff changeset
9 Free Software Foundation; either version 3 of the License, or (at your
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 5823
diff changeset
10 option) any later version.
2928
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
11
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
12 Octave is distributed in the hope that it will be useful, but WITHOUT
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
13 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
14 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
15 for more details.
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
16
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
17 You should have received a copy of the GNU General Public License
7016
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 5823
diff changeset
18 along with Octave; see the file COPYING. If not, see
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 5823
diff changeset
19 <http://www.gnu.org/licenses/>.
2928
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
20
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
21 */
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
22
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
23 #ifdef HAVE_CONFIG_H
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
24 #include <config.h>
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
25 #endif
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
26
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
27 #include "CmplxLU.h"
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
28 #include "dbleLU.h"
7789
82be108cc558 First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents: 7515
diff changeset
29 #include "fCmplxLU.h"
82be108cc558 First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents: 7515
diff changeset
30 #include "floatLU.h"
7515
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
31 #include "SparseCmplxLU.h"
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
32 #include "SparsedbleLU.h"
2928
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
33
15039
e753177cde93 maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents: 14854
diff changeset
34 #include "defun.h"
2928
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
35 #include "error.h"
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
36 #include "gripes.h"
20940
48b2ad5ee801 maint: Rename oct-obj.[cc|h] to ovl.[cc|h] for clarity.
Rik <rik@octave.org>
parents: 20921
diff changeset
37 #include "ovl.h"
2928
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
38 #include "utils.h"
7515
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
39 #include "ov-re-sparse.h"
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
40 #include "ov-cx-sparse.h"
2928
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
41
8869
c3b743b1b1c6 preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents: 8367
diff changeset
42 template <class MT>
c3b743b1b1c6 preset triangular type if possible for lu and qr outputs
Jaroslav Hajek <highegg@gmail.com>
parents: 8367
diff changeset
43 static octave_value
9715
9f27172fbd1e auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents: 9709
diff changeset
44 get_lu_l (const base_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 ();
9f27172fbd1e auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents: 9709
diff changeset
47 if (L.is_square ())
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
9715
9f27172fbd1e auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents: 9709
diff changeset
53 template <class MT>
9f27172fbd1e auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents: 9709
diff changeset
54 static octave_value
9f27172fbd1e auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents: 9709
diff changeset
55 get_lu_u (const base_lu<MT>& fact)
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 ();
9f27172fbd1e auto-set MatrixType from certain functions
Jaroslav Hajek <highegg@gmail.com>
parents: 9709
diff changeset
58 if (U.is_square () && fact.regular ())
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,
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
65 "-*- texinfo -*-\n\
20853
1142cf6abc0d 2015 Code Sprint: remove class of function from docstring for all C++ files.
Rik <rik@octave.org>
parents: 20816
diff changeset
66 @deftypefn {} {[@var{L}, @var{U}] =} lu (@var{A})\n\
1142cf6abc0d 2015 Code Sprint: remove class of function from docstring for all C++ files.
Rik <rik@octave.org>
parents: 20816
diff changeset
67 @deftypefnx {} {[@var{L}, @var{U}, @var{P}] =} lu (@var{A})\n\
1142cf6abc0d 2015 Code Sprint: remove class of function from docstring for all C++ files.
Rik <rik@octave.org>
parents: 20816
diff changeset
68 @deftypefnx {} {[@var{L}, @var{U}, @var{P}, @var{Q}] =} lu (@var{S})\n\
1142cf6abc0d 2015 Code Sprint: remove class of function from docstring for all C++ files.
Rik <rik@octave.org>
parents: 20816
diff changeset
69 @deftypefnx {} {[@var{L}, @var{U}, @var{P}, @var{Q}, @var{R}] =} lu (@var{S})\n\
1142cf6abc0d 2015 Code Sprint: remove class of function from docstring for all C++ files.
Rik <rik@octave.org>
parents: 20816
diff changeset
70 @deftypefnx {} {[@dots{}] =} lu (@var{S}, @var{thres})\n\
1142cf6abc0d 2015 Code Sprint: remove class of function from docstring for all C++ files.
Rik <rik@octave.org>
parents: 20816
diff changeset
71 @deftypefnx {} {@var{y} =} lu (@dots{})\n\
1142cf6abc0d 2015 Code Sprint: remove class of function from docstring for all C++ files.
Rik <rik@octave.org>
parents: 20816
diff changeset
72 @deftypefnx {} {[@dots{}] =} lu (@dots{}, \"vector\")\n\
3372
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
73 @cindex LU decomposition\n\
20172
4f45eaf83908 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19864
diff changeset
74 Compute the LU@tie{}decomposition of @var{A}.\n\
4f45eaf83908 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19864
diff changeset
75 \n\
4f45eaf83908 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19864
diff changeset
76 If @var{A} is full subroutines from @sc{lapack} are used and if @var{A} is\n\
4f45eaf83908 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19864
diff changeset
77 sparse then @sc{umfpack} is used.\n\
4f45eaf83908 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19864
diff changeset
78 \n\
4f45eaf83908 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19864
diff changeset
79 The result is returned in a permuted form, according to the optional return\n\
11553
01f703952eff Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
80 value @var{P}. For example, given the matrix @code{a = [1, 2; 3, 4]},\n\
3372
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
81 \n\
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
82 @example\n\
11553
01f703952eff Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
83 [l, u, p] = lu (@var{a})\n\
3372
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
84 @end example\n\
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
85 \n\
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
86 @noindent\n\
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
87 returns\n\
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
88 \n\
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
89 @example\n\
9064
7c02ec148a3c Check grammar on all .cc files
Rik <rdrider0-list@yahoo.com>
parents: 8969
diff changeset
90 @group\n\
3372
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
91 l =\n\
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
92 \n\
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
93 1.00000 0.00000\n\
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
94 0.33333 1.00000\n\
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
95 \n\
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
96 u =\n\
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
97 \n\
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
98 3.00000 4.00000\n\
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
99 0.00000 0.66667\n\
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
100 \n\
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
101 p =\n\
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
102 \n\
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
103 0 1\n\
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
104 1 0\n\
9064
7c02ec148a3c Check grammar on all .cc files
Rik <rdrider0-list@yahoo.com>
parents: 8969
diff changeset
105 @end group\n\
3372
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
106 @end example\n\
4329
d53c33d93440 [project @ 2003-02-18 20:00:48 by jwe]
jwe
parents: 3587
diff changeset
107 \n\
5457
c6dc1ccd83a9 [project @ 2005-09-19 19:23:35 by jwe]
jwe
parents: 5307
diff changeset
108 The matrix is not required to be square.\n\
7515
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
109 \n\
11553
01f703952eff Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
110 When called with two or three output arguments and a spare input matrix,\n\
01f703952eff Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
111 @code{lu} does not attempt to perform sparsity preserving column\n\
9064
7c02ec148a3c Check grammar on all .cc files
Rik <rdrider0-list@yahoo.com>
parents: 8969
diff changeset
112 permutations. Called with a fourth output argument, the sparsity\n\
7515
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
113 preserving column transformation @var{Q} is returned, such that\n\
11553
01f703952eff Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
114 @code{@var{P} * @var{A} * @var{Q} = @var{L} * @var{U}}.\n\
7515
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
115 \n\
11553
01f703952eff Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
116 Called with a fifth output argument and a sparse input matrix,\n\
01f703952eff Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
117 @code{lu} attempts to use a scaling factor @var{R} on the input matrix\n\
10846
a4f482e66b65 Grammarcheck more of the documentation.
Rik <octave@nomad.inbox5.com>
parents: 10840
diff changeset
118 such that\n\
11553
01f703952eff Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
119 @code{@var{P} * (@var{R} \\ @var{A}) * @var{Q} = @var{L} * @var{U}}.\n\
9065
8207b833557f Cleanup documentation for arith.texi, linalg.texi, nonlin.texi
Rik <rdrider0-list@yahoo.com>
parents: 9064
diff changeset
120 This typically leads to a sparser and more stable factorization.\n\
7515
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
121 \n\
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
122 An additional input argument @var{thres}, that defines the pivoting\n\
9064
7c02ec148a3c Check grammar on all .cc files
Rik <rdrider0-list@yahoo.com>
parents: 8969
diff changeset
123 threshold can be given. @var{thres} can be a scalar, in which case\n\
11553
01f703952eff Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
124 it defines the @sc{umfpack} pivoting tolerance for both symmetric and\n\
01f703952eff Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
125 unsymmetric cases. If @var{thres} is a 2-element vector, then the first\n\
01f703952eff Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
126 element defines the pivoting tolerance for the unsymmetric @sc{umfpack}\n\
01f703952eff Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
127 pivoting strategy and the second for the symmetric strategy. By default,\n\
01f703952eff Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
128 the values defined by @code{spparms} are used ([0.1, 0.001]).\n\
7515
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
129 \n\
17281
bc924baa2c4e doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents: 17268
diff changeset
130 Given the string argument @qcode{\"vector\"}, @code{lu} returns the values\n\
bc924baa2c4e doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents: 17268
diff changeset
131 of @var{P} and @var{Q} as vector values, such that for full matrix,\n\
bc924baa2c4e doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents: 17268
diff changeset
132 @code{@var{A} (@var{P},:) = @var{L} * @var{U}}, and @code{@var{R}(@var{P},:)\n\
bc924baa2c4e doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents: 17268
diff changeset
133 * @var{A} (:, @var{Q}) = @var{L} * @var{U}}.\n\
7515
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
134 \n\
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
135 With two output arguments, returns the permuted forms of the upper and\n\
11553
01f703952eff Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
136 lower triangular matrices, such that @code{@var{A} = @var{L} * @var{U}}.\n\
9209
923c7cb7f13f Simplify TeXinfo files by eliminating redundant @iftex followed by @tex construction.
Rik <rdrider0-list@yahoo.com>
parents: 9065
diff changeset
137 With one output argument @var{y}, then the matrix returned by the @sc{lapack}\n\
11553
01f703952eff Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
138 routines is returned. If the input matrix is sparse then the matrix @var{L}\n\
01f703952eff Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
139 is embedded into @var{U} to give a return value similar to the full case.\n\
01f703952eff Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
140 For both full and sparse matrices, @code{lu} loses the permutation\n\
7515
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
141 information.\n\
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 18849
diff changeset
142 @seealso{luupdate, ilu, chol, hess, qr, qz, schur, svd}\n\
3372
f16c2ce14886 [project @ 1999-11-23 19:07:09 by jwe]
jwe
parents: 3014
diff changeset
143 @end deftypefn")
2928
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
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 ();
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
146 bool issparse = (nargin > 0 && args(0).is_sparse_type ());
2928
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
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
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
150
7515
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
151 bool vecout = false;
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
152 Matrix thres;
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
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
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
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
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
168 if (! issparse)
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
169 error ("lu: can not define pivoting threshold THRES for full matrices");
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10094
diff changeset
170 Matrix tmp = args(n++).matrix_value ();
7515
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
171
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
172 if (tmp.numel () == 1)
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10094
diff changeset
173 {
20555
f90c8372b7ba eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents: 20228
diff changeset
174 thres.resize (1,2);
f90c8372b7ba eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents: 20228
diff changeset
175 thres(0) = tmp(0);
f90c8372b7ba eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents: 20228
diff changeset
176 thres(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)
f90c8372b7ba eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents: 20228
diff changeset
179 thres = tmp;
f90c8372b7ba eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents: 20228
diff changeset
180 else
20711
7b608fadc663 Make error messages more specific about the variable and problem encountered.
Rik <rik@octave.org>
parents: 20555
diff changeset
181 error ("lu: THRES 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
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
185 octave_value_list retval;
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
186 bool scale = (nargout == 5);
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
187
2928
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
188 octave_value arg = args(0);
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
189
5275
23b37da9fd5b [project @ 2005-04-08 16:07:35 by jwe]
jwe
parents: 4329
diff changeset
190 octave_idx_type nr = arg.rows ();
23b37da9fd5b [project @ 2005-04-08 16:07:35 by jwe]
jwe
parents: 4329
diff changeset
191 octave_idx_type nc = arg.columns ();
2928
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
192
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
193 int arg_is_empty = empty_arg ("lu", nr, nc);
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
194
7515
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
195 if (issparse)
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
196 {
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
197 if (arg_is_empty < 0)
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
198 return octave_value_list ();
7515
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
199 else if (arg_is_empty > 0)
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10094
diff changeset
200 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
201
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
202 if (arg.is_real_type ())
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10094
diff changeset
203 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10094
diff changeset
204 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
205
18678
6113e0c6920b maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents: 18497
diff changeset
206 if (nargout < 4)
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10094
diff changeset
207 {
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
208 ColumnVector Qinit (nc);
18425
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
209 for (octave_idx_type i = 0; i < nc; i++)
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
210 Qinit(i) = i;
18425
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
211 SparseLU fact (m, Qinit, thres, false, 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)
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
219 = octave_value (
19593
446c46af4b42 strip trailing whitespace from most source files
John W. Eaton <jwe@octave.org>
parents: 18488
diff changeset
220 fact.U () * fact.Pc_mat ().transpose (),
18425
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
221 MatrixType (MatrixType::Permuted_Upper,
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
222 nc, fact.col_perm ()));
7515
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
223
18425
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
224 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
225 SparseMatrix L = fact.L ();
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
226
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
227 if (nargout == 2)
18425
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
228 retval(0)
18678
6113e0c6920b maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents: 18497
diff changeset
229 = octave_value (P.transpose () * L,
18425
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
230 MatrixType (MatrixType::Permuted_Lower,
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
231 nr, fact.row_perm ()));
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
232 else
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
233 {
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
234 retval(0) = L;
18678
6113e0c6920b maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents: 18497
diff changeset
235 if (vecout)
18425
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
236 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
237 else
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
238 retval(2) = P;
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 }
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
242 else
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
243 {
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
244 retval.resize (scale ? 5 : 4);
19861
19755f4fc851 maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19743
diff changeset
245 SparseLU fact (m, thres, scale);
2928
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
246
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
247 retval(0) = octave_value (fact.L (),
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
248 MatrixType (MatrixType::Lower));
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
249 retval(1) = octave_value (fact.U (),
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
250 MatrixType (MatrixType::Upper));
7515
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
251
19861
19755f4fc851 maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19743
diff changeset
252 if (vecout)
19755f4fc851 maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19743
diff changeset
253 {
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
254 retval(2) = fact.Pr_vec ();
19861
19755f4fc851 maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19743
diff changeset
255 retval(3) = fact.Pc_vec ();
19755f4fc851 maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19743
diff changeset
256 }
19755f4fc851 maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19743
diff changeset
257 else
19755f4fc851 maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19743
diff changeset
258 {
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
259 retval(2) = fact.Pr_mat ();
19861
19755f4fc851 maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19743
diff changeset
260 retval(3) = fact.Pc_mat ();
19755f4fc851 maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19743
diff changeset
261 }
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
262
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
263 if (scale)
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
264 retval(4) = fact.R ();
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10094
diff changeset
265 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10094
diff changeset
266 }
7515
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
267 else if (arg.is_complex_type ())
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10094
diff changeset
268 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10094
diff changeset
269 SparseComplexMatrix m = arg.sparse_complex_matrix_value ();
2928
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
270
18678
6113e0c6920b maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents: 18497
diff changeset
271 if (nargout < 4)
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10094
diff changeset
272 {
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
273 ColumnVector Qinit (nc);
18425
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
274 for (octave_idx_type i = 0; i < nc; i++)
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
275 Qinit(i) = i;
18425
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
276 SparseComplexLU fact (m, Qinit, thres, false, true);
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
277
18678
6113e0c6920b maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents: 18497
diff changeset
278 if (nargout < 2)
19861
19755f4fc851 maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19743
diff changeset
279 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
280 else
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
281 {
20921
4d3daf7e43f3 eliminate trailing whitespace in source files
John W. Eaton <jwe@octave.org>
parents: 20909
diff changeset
282 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
283 retval(1)
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
284 = octave_value (
19593
446c46af4b42 strip trailing whitespace from most source files
John W. Eaton <jwe@octave.org>
parents: 18488
diff changeset
285 fact.U () * fact.Pc_mat ().transpose (),
18425
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
286 MatrixType (MatrixType::Permuted_Upper,
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
287 nc, fact.col_perm ()));
7515
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
288
18425
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
289 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
290 SparseComplexMatrix L = fact.L ();
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
291 if (nargout == 2)
19861
19755f4fc851 maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19743
diff changeset
292 retval(0)
19755f4fc851 maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19743
diff changeset
293 = octave_value (P.transpose () * L,
19755f4fc851 maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19743
diff changeset
294 MatrixType (MatrixType::Permuted_Lower,
19755f4fc851 maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19743
diff changeset
295 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
296 else
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
297 {
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
298 retval(0) = L;
18678
6113e0c6920b maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents: 18497
diff changeset
299 if (vecout)
18425
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
300 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
301 else
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
302 retval(2) = P;
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
303 }
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
304 }
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 else
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
307 {
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
308 retval.resize (scale ? 5 : 4);
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
309 SparseComplexLU fact (m, thres, scale);
7515
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
310
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
311 retval(0) = octave_value (fact.L (),
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
312 MatrixType (MatrixType::Lower));
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
313 retval(1) = octave_value (fact.U (),
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
314 MatrixType (MatrixType::Upper));
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
315
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
316 if (vecout)
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
317 {
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
318 retval(2) = fact.Pr_vec ();
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
319 retval(3) = fact.Pc_vec ();
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
320 }
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
321 else
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
322 {
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
323 retval(2) = fact.Pr_mat ();
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
324 retval(3) = fact.Pc_mat ();
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
325 }
7515
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
326
19861
19755f4fc851 maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19743
diff changeset
327 if (scale)
19755f4fc851 maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19743
diff changeset
328 retval(4) = fact.R ();
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10094
diff changeset
329 }
18427
93c019f00e59 maint: whitespace fixes for libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 18426
diff changeset
330
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10094
diff changeset
331 }
7515
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
332 else
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10094
diff changeset
333 gripe_wrong_type_arg ("lu", arg);
2928
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
334 }
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
335 else
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
336 {
7515
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
337 if (arg_is_empty < 0)
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
338 return octave_value_list ();
7515
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
339 else if (arg_is_empty > 0)
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
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
342 if (arg.is_real_type ())
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
20555
f90c8372b7ba eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents: 20228
diff changeset
348 FloatLU 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
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
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
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
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
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
369 retval = ovl (get_lu_l (fact), get_lu_u (fact),
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
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
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
372 retval = ovl (get_lu_l (fact), get_lu_u (fact),
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
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
20555
f90c8372b7ba eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents: 20228
diff changeset
382 LU 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
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
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
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
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
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
403 retval = ovl (get_lu_l (fact), get_lu_u (fact),
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
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
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
406 retval = ovl (get_lu_l (fact), get_lu_u (fact),
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
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 }
7515
f3c00dc0912b Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents: 7017
diff changeset
413 else if (arg.is_complex_type ())
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
20555
f90c8372b7ba eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents: 20228
diff changeset
419 FloatComplexLU 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
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
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
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
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
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
440 retval = ovl (get_lu_l (fact), get_lu_u (fact),
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
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
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
443 retval = ovl (get_lu_l (fact), get_lu_u (fact),
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
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
20555
f90c8372b7ba eliminate many more simple uses of error_state
John W. Eaton <jwe@octave.org>
parents: 20228
diff changeset
453 ComplexLU 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
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
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
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
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
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
474 retval = ovl (get_lu_l (fact), get_lu_u (fact),
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
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
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
477 retval = ovl (get_lu_l (fact), get_lu_u (fact),
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
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
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10094
diff changeset
485 gripe_wrong_type_arg ("lu", arg);
2928
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
486 }
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
487
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
488 return retval;
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
489 }
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
490
295f037b4b3e [project @ 1997-05-05 05:32:33 by jwe]
jwe
parents:
diff changeset
491 /*
7814
87865ed7405f Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents: 7789
diff changeset
492 %!assert(lu ([1, 2; 3, 4]), [3, 4; 1/3, 2/3], eps);
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
14501
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14361
diff changeset
542 %!error lu ()
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14361
diff changeset
543 %!error <can not define pivoting threshold> lu ([1, 2; 3, 4], 2)
18425
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
544
18488
4c2465444a96 Use %!testif HAVE_UMFPACK for sparse lu test added in cset 2a45b6b87bee
Mike Miller <mtmiller@ieee.org>
parents: 18427
diff changeset
545 %!testif HAVE_UMFPACK
18425
2a45b6b87bee correct numerical errors in sparse LU factorization (bug #41116).
Michael C. Grant <mcg@cvxr.com>
parents: 17787
diff changeset
546 %! 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
547 %! 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
548 %! 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
549 %! B = sparse (Bi, Bj, Bv);
93c019f00e59 maint: whitespace fixes for libinterp/corefcn/lu.cc
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 18426
diff changeset
550 %! [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
551 %! 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
552 %! [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
553 %! 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
554 %! [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
555 %! 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
556
14501
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14361
diff changeset
557 */
7814
87865ed7405f Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents: 7789
diff changeset
558
9708
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
559 static
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
560 bool check_lu_dims (const octave_value& l, const octave_value& u,
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
561 const octave_value& p)
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
562 {
18100
6a71e5030df5 Follow coding convention of defining and initializing only 1 variable per line in liboctinterp.
Rik <rik@octave.org>
parents: 17787
diff changeset
563 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
564 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
565 octave_idx_type n = u.columns ();
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
566
9708
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
567 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
568 && 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
569 && (p.is_undefined () || p.rows () == m));
9708
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
570 }
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
571
15039
e753177cde93 maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents: 14854
diff changeset
572 DEFUN (luupdate, args, ,
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
573 "-*- texinfo -*-\n\
20853
1142cf6abc0d 2015 Code Sprint: remove class of function from docstring for all C++ files.
Rik <rik@octave.org>
parents: 20816
diff changeset
574 @deftypefn {} {[@var{L}, @var{U}] =} luupdate (@var{L}, @var{U}, @var{x}, @var{y})\n\
1142cf6abc0d 2015 Code Sprint: remove class of function from docstring for all C++ files.
Rik <rik@octave.org>
parents: 20816
diff changeset
575 @deftypefnx {} {[@var{L}, @var{U}, @var{P}] =} luupdate (@var{L}, @var{U}, @var{P}, @var{x}, @var{y})\n\
9709
f8e2e9fdaa8f document luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9708
diff changeset
576 Given an LU@tie{}factorization of a real or complex matrix\n\
f8e2e9fdaa8f document luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9708
diff changeset
577 @w{@var{A} = @var{L}*@var{U}}, @var{L}@tie{}lower unit trapezoidal and\n\
f8e2e9fdaa8f document luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9708
diff changeset
578 @var{U}@tie{}upper trapezoidal, return the LU@tie{}factorization\n\
f8e2e9fdaa8f document luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9708
diff changeset
579 of @w{@var{A} + @var{x}*@var{y}.'}, where @var{x} and @var{y} are\n\
f8e2e9fdaa8f document luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9708
diff changeset
580 column vectors (rank-1 update) or matrices with equal number of columns\n\
f8e2e9fdaa8f document luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9708
diff changeset
581 (rank-k update).\n\
20172
4f45eaf83908 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19864
diff changeset
582 \n\
4f45eaf83908 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19864
diff changeset
583 Optionally, row-pivoted updating can be used by supplying a row permutation\n\
4f45eaf83908 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19864
diff changeset
584 (pivoting) matrix @var{P}; in that case, an updated permutation matrix is\n\
4f45eaf83908 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19864
diff changeset
585 returned. Note that if @var{L}, @var{U}, @var{P} is a pivoted\n\
4f45eaf83908 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19864
diff changeset
586 LU@tie{}factorization as obtained by @code{lu}:\n\
9709
f8e2e9fdaa8f document luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9708
diff changeset
587 \n\
f8e2e9fdaa8f document luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9708
diff changeset
588 @example\n\
14360
97883071e8e4 doc: Correct off-by-1 spacings in all .cc docstrings
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
589 [@var{L}, @var{U}, @var{P}] = lu (@var{A});\n\
9709
f8e2e9fdaa8f document luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9708
diff changeset
590 @end example\n\
f8e2e9fdaa8f document luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9708
diff changeset
591 \n\
10846
a4f482e66b65 Grammarcheck more of the documentation.
Rik <octave@nomad.inbox5.com>
parents: 10840
diff changeset
592 @noindent\n\
17268
1c21f264d26f doc: Rename @xcode macro to @tcode (transpose code) for clarity.
Rik <rik@octave.org>
parents: 16920
diff changeset
593 then a factorization of @tcode{@var{A}+@var{x}*@var{y}.'} can be obtained\n\
11553
01f703952eff Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
594 either as\n\
9709
f8e2e9fdaa8f document luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9708
diff changeset
595 \n\
f8e2e9fdaa8f document luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9708
diff changeset
596 @example\n\
14360
97883071e8e4 doc: Correct off-by-1 spacings in all .cc docstrings
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
597 [@var{L1}, @var{U1}] = lu (@var{L}, @var{U}, @var{P}*@var{x}, @var{y})\n\
9709
f8e2e9fdaa8f document luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9708
diff changeset
598 @end example\n\
f8e2e9fdaa8f document luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9708
diff changeset
599 \n\
10846
a4f482e66b65 Grammarcheck more of the documentation.
Rik <octave@nomad.inbox5.com>
parents: 10840
diff changeset
600 @noindent\n\
9709
f8e2e9fdaa8f document luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9708
diff changeset
601 or\n\
f8e2e9fdaa8f document luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9708
diff changeset
602 \n\
f8e2e9fdaa8f document luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9708
diff changeset
603 @example\n\
14360
97883071e8e4 doc: Correct off-by-1 spacings in all .cc docstrings
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
604 [@var{L1}, @var{U1}, @var{P1}] = lu (@var{L}, @var{U}, @var{P}, @var{x}, @var{y})\n\
9709
f8e2e9fdaa8f document luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9708
diff changeset
605 @end example\n\
f8e2e9fdaa8f document luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9708
diff changeset
606 \n\
11553
01f703952eff Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
607 The first form uses the unpivoted algorithm, which is faster, but less\n\
01f703952eff Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
608 stable. The second form uses a slower pivoted algorithm, which is more\n\
01f703952eff Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
609 stable.\n\
9709
f8e2e9fdaa8f document luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9708
diff changeset
610 \n\
20172
4f45eaf83908 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19864
diff changeset
611 The matrix case is done as a sequence of rank-1 updates; thus, for large\n\
4f45eaf83908 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19864
diff changeset
612 enough k, it will be both faster and more accurate to recompute the\n\
4f45eaf83908 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19864
diff changeset
613 factorization from scratch.\n\
16920
53eaa83e4181 doc: Add seealso links between various factorization forms.
Rik <rik@octave.org>
parents: 15195
diff changeset
614 @seealso{lu, cholupdate, qrupdate}\n\
9724
f22bbc5d56e9 Fix various incorrect usages of TeXinfo deffn and deftypefn macros
Rik <rdrider0-list@yahoo.com>
parents: 9715
diff changeset
615 @end deftypefn")
9708
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
616 {
20816
b16bcd7a2a33 Use int rather than octave_idx_type for nargin data type.
Rik <rik@octave.org>
parents: 20812
diff changeset
617 int nargin = args.length ();
9708
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
618
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
619 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
620 print_usage ();
9708
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
621
20812
d9ca869ca124 maint: Clean-up more instances of print_usage().
Rik <rik@octave.org>
parents: 20802
diff changeset
622 bool pivoted = (nargin == 5);
d9ca869ca124 maint: Clean-up more instances of print_usage().
Rik <rik@octave.org>
parents: 20802
diff changeset
623
9708
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
624 octave_value argl = args(0);
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
625 octave_value argu = args(1);
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
626 octave_value argp = pivoted ? args(2) : octave_value ();
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
627 octave_value argx = args(2 + pivoted);
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
628 octave_value argy = args(3 + pivoted);
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
629
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
630 if (! (argl.is_numeric_type () && argu.is_numeric_type ()
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
631 && argx.is_numeric_type () && argy.is_numeric_type ()
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
632 && (! pivoted || argp.is_perm_matrix ())))
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
633 error ("luupdate: L, U, X, and Y must be numeric");
9708
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
634
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
635 if (! check_lu_dims (argl, argu, argp))
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
636 error ("luupdate: dimension mismatch");
9708
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
637
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
638 PermMatrix P = (pivoted
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
639 ? argp.perm_matrix_value ()
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
640 : PermMatrix::eye (argl.rows ()));
9708
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
641
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
642 if (argl.is_real_type () && argu.is_real_type ()
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
643 && argx.is_real_type () && argy.is_real_type ())
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
644 {
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
645 // all real case
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
646 if (argl.is_single_type () || argu.is_single_type ()
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
647 || argx.is_single_type () || argy.is_single_type ())
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
648 {
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
649 FloatMatrix L = argl.float_matrix_value ();
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
650 FloatMatrix U = argu.float_matrix_value ();
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
651 FloatMatrix x = argx.float_matrix_value ();
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
652 FloatMatrix y = argy.float_matrix_value ();
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11572
diff changeset
653
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
654 FloatLU fact (L, U, P);
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
655 if (pivoted)
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
656 fact.update_piv (x, y);
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
657 else
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
658 fact.update (x, y);
9708
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
659
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
660 if (pivoted)
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20892
diff changeset
661 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
662 else
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20892
diff changeset
663 return ovl (get_lu_l (fact), get_lu_u (fact));
9708
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
664 }
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
665 else
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
666 {
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
667 Matrix L = argl.matrix_value ();
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
668 Matrix U = argu.matrix_value ();
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
669 Matrix x = argx.matrix_value ();
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
670 Matrix y = argy.matrix_value ();
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
671
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
672 LU fact (L, U, P);
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
673 if (pivoted)
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
674 fact.update_piv (x, y);
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
675 else
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
676 fact.update (x, y);
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
677
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
678 if (pivoted)
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20892
diff changeset
679 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
680 else
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20892
diff changeset
681 return ovl (get_lu_l (fact), get_lu_u (fact));
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
682 }
9708
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
683 }
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
684 else
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
685 {
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
686 // complex case
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
687 if (argl.is_single_type () || argu.is_single_type ()
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
688 || argx.is_single_type () || argy.is_single_type ())
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
689 {
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
690 FloatComplexMatrix L = argl.float_complex_matrix_value ();
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
691 FloatComplexMatrix U = argu.float_complex_matrix_value ();
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
692 FloatComplexMatrix x = argx.float_complex_matrix_value ();
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
693 FloatComplexMatrix y = argy.float_complex_matrix_value ();
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
694
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
695 FloatComplexLU fact (L, U, P);
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
696 if (pivoted)
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
697 fact.update_piv (x, y);
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
698 else
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
699 fact.update (x, y);
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
700
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
701 if (pivoted)
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20892
diff changeset
702 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
703 else
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20892
diff changeset
704 return ovl (get_lu_l (fact), get_lu_u (fact));
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
705 }
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
706 else
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
707 {
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
708 ComplexMatrix L = argl.complex_matrix_value ();
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
709 ComplexMatrix U = argu.complex_matrix_value ();
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
710 ComplexMatrix x = argx.complex_matrix_value ();
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
711 ComplexMatrix y = argy.complex_matrix_value ();
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
712
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
713 ComplexLU fact (L, U, P);
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
714 if (pivoted)
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
715 fact.update_piv (x, y);
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
716 else
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
717 fact.update (x, y);
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
718
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
719 if (pivoted)
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20892
diff changeset
720 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
721 else
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20892
diff changeset
722 return ovl (get_lu_l (fact), get_lu_u (fact));
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
723 }
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
724 }
9708
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
725 }
6f3ffe11d926 implement luupdate
Jaroslav Hajek <highegg@gmail.com>
parents: 9209
diff changeset
726
7814
87865ed7405f Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents: 7789
diff changeset
727 /*
10080
cf70ee43077c add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents: 9724
diff changeset
728 %!shared A, u, v, Ac, uc, vc
cf70ee43077c add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents: 9724
diff changeset
729 %! A = [0.091364 0.613038 0.999083;
cf70ee43077c add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents: 9724
diff changeset
730 %! 0.594638 0.425302 0.603537;
cf70ee43077c add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents: 9724
diff changeset
731 %! 0.383594 0.291238 0.085574;
cf70ee43077c add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents: 9724
diff changeset
732 %! 0.265712 0.268003 0.238409;
cf70ee43077c add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents: 9724
diff changeset
733 %! 0.669966 0.743851 0.445057 ];
cf70ee43077c add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents: 9724
diff changeset
734 %!
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11572
diff changeset
735 %! u = [0.85082;
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11572
diff changeset
736 %! 0.76426;
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11572
diff changeset
737 %! 0.42883;
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11572
diff changeset
738 %! 0.53010;
10080
cf70ee43077c add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents: 9724
diff changeset
739 %! 0.80683 ];
cf70ee43077c add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents: 9724
diff changeset
740 %!
cf70ee43077c add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents: 9724
diff changeset
741 %! v = [0.98810;
cf70ee43077c add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents: 9724
diff changeset
742 %! 0.24295;
cf70ee43077c add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents: 9724
diff changeset
743 %! 0.43167 ];
cf70ee43077c add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents: 9724
diff changeset
744 %!
cf70ee43077c add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents: 9724
diff changeset
745 %! 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
746 %! 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
747 %! 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
748 %! 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
749 %! 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
750 %!
cf70ee43077c add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents: 9724
diff changeset
751 %! 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
752 %! 0.13141 + 0.43708i;
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14361
diff changeset
753 %! 0.29808 + 0.08789i;
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14361
diff changeset
754 %! 0.69821 + 0.38844i;
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14361
diff changeset
755 %! 0.74871 + 0.25821i ];
10080
cf70ee43077c add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents: 9724
diff changeset
756 %!
cf70ee43077c add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents: 9724
diff changeset
757 %! 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
758 %! 0.20820 + 0.93090i;
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14361
diff changeset
759 %! 0.86184 + 0.34689i ];
10080
cf70ee43077c add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents: 9724
diff changeset
760 %!
cf70ee43077c add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents: 9724
diff changeset
761
10094
ab1101011a6d lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents: 10080
diff changeset
762 %!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
763 %! [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
764 %! [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
765 %! 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
766 %! 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
767 %! 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
768 %!
10094
ab1101011a6d lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents: 10080
diff changeset
769 %!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
770 %! [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
771 %! [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
772 %! 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
773 %! 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
774 %! 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
775
10094
ab1101011a6d lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents: 10080
diff changeset
776 %!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
777 %! [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
778 %! [L,U] = 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
779 %! 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
780 %! 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
781 %! assert (norm (vec (P'*L*U - single (A) - single (u)*single (v).'), Inf) < norm (single (A))*1e1*eps ("single"));
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11572
diff changeset
782 %!
10094
ab1101011a6d lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents: 10080
diff changeset
783 %!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
784 %! [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
785 %! [L,U] = 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
786 %! 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
787 %! 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
788 %! assert (norm (vec (P'*L*U - single (Ac) - single (uc)*single (vc).'), Inf) < norm (single (Ac))*1e1*eps ("single"));
10080
cf70ee43077c add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents: 9724
diff changeset
789
10094
ab1101011a6d lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents: 10080
diff changeset
790 %!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
791 %! [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
792 %! [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
793 %! 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
794 %! 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
795 %! 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
796 %!
10094
ab1101011a6d lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents: 10080
diff changeset
797 %!testif HAVE_QRUPDATE_LUU
18849
aa9ca67f09fb make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents: 18678
diff changeset
798 %! [L,U,P] = lu (A);
aa9ca67f09fb make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents: 18678
diff changeset
799 %! [~,ordcols] = max (P,[],1);
aa9ca67f09fb make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents: 18678
diff changeset
800 %! [~,ordrows] = max (P,[],2);
aa9ca67f09fb make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents: 18678
diff changeset
801 %! P1 = eye (size(P))(:,ordcols);
aa9ca67f09fb make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents: 18678
diff changeset
802 %! P2 = eye (size(P))(ordrows,:);
aa9ca67f09fb make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents: 18678
diff changeset
803 %! assert(P1 == P);
aa9ca67f09fb make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents: 18678
diff changeset
804 %! assert(P2 == P);
aa9ca67f09fb make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents: 18678
diff changeset
805 %! [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
806 %! [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
807 %! [L,U,P2] = luupdate (L,U,P2,u,v);
aa9ca67f09fb make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents: 18678
diff changeset
808 %! assert(P1 == P);
aa9ca67f09fb make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents: 18678
diff changeset
809 %! assert(P2 == P);
aa9ca67f09fb make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents: 18678
diff changeset
810 %!
aa9ca67f09fb make all permutation matrices column permutations (bug #42418)
David Spies <dnspies@gmail.com>
parents: 18678
diff changeset
811 %!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
812 %! [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
813 %! [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
814 %! 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
815 %! 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
816 %! 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
817
10094
ab1101011a6d lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents: 10080
diff changeset
818 %!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
819 %! [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
820 %! [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
821 %! 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
822 %! 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
823 %! assert (norm (vec (P'*L*U - single (A) - single (u)*single (v).'), Inf) < norm (single (A))*1e1*eps ("single"));
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11572
diff changeset
824 %!
10094
ab1101011a6d lu.cc: avoid test failures if HAVE_QRUPDATE_LUU is not defined
John W. Eaton <jwe@octave.org>
parents: 10080
diff changeset
825 %!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
826 %! [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
827 %! [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
828 %! 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
829 %! 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
830 %! assert (norm (vec (P'*L*U - single (Ac) - single (uc)*single (vc).'), Inf) < norm (single (Ac))*1e1*eps ("single"));
10080
cf70ee43077c add tests for LU updating
Jaroslav Hajek <highegg@gmail.com>
parents: 9724
diff changeset
831 */