GNU Octave  9.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
lu.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 //
3 // Copyright (C) 1996-2024 The Octave Project Developers
4 //
5 // See the file COPYRIGHT.md in the top-level directory of this
6 // distribution or <https://octave.org/copyright/>.
7 //
8 // This file is part of Octave.
9 //
10 // Octave is free software: you can redistribute it and/or modify it
11 // under the terms of the GNU General Public License as published by
12 // the Free Software Foundation, either version 3 of the License, or
13 // (at your option) any later version.
14 //
15 // Octave is distributed in the hope that it will be useful, but
16 // WITHOUT ANY WARRANTY; without even the implied warranty of
17 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 // GNU General Public License for more details.
19 //
20 // You should have received a copy of the GNU General Public License
21 // along with Octave; see the file COPYING. If not, see
22 // <https://www.gnu.org/licenses/>.
23 //
24 ////////////////////////////////////////////////////////////////////////
25 
26 #if defined (HAVE_CONFIG_H)
27 # include "config.h"
28 #endif
29 
30 #include "lu.h"
31 #include "sparse-lu.h"
32 
33 #include "defun.h"
34 #include "error.h"
35 #include "errwarn.h"
36 #include "ovl.h"
37 #include "ov-re-sparse.h"
38 #include "ov-cx-sparse.h"
39 
41 
42 template <typename MT>
43 static octave_value
44 get_lu_l (const math::lu<MT>& fact)
45 {
46  MT L = fact.L ();
47  if (L.issquare ())
49  else
50  return L;
51 }
52 
53 template <typename MT>
54 static octave_value
55 get_lu_u (const math::lu<MT>& fact)
56 {
57  MT U = fact.U ();
58  if (U.issquare () && fact.regular ())
60  else
61  return U;
62 }
63 
64 DEFUN (lu, args, nargout,
65  doc: /* -*- texinfo -*-
66 @deftypefn {} {[@var{L}, @var{U}] =} lu (@var{A})
67 @deftypefnx {} {[@var{L}, @var{U}, @var{P}] =} lu (@var{A})
68 @deftypefnx {} {[@var{L}, @var{U}, @var{P}, @var{Q}] =} lu (@var{S})
69 @deftypefnx {} {[@var{L}, @var{U}, @var{P}, @var{Q}, @var{R}] =} lu (@var{S})
70 @deftypefnx {} {[@dots{}] =} lu (@var{S}, @var{thresh})
71 @deftypefnx {} {@var{y} =} lu (@dots{})
72 @deftypefnx {} {[@dots{}] =} lu (@dots{}, "vector")
73 @cindex LU decomposition
74 Compute the LU@tie{}decomposition of @var{A}.
75 
76 If @var{A} is full then subroutines from @sc{lapack} are used, and if
77 @var{A} is sparse then @sc{umfpack} is used.
78 
79 The result is returned in a permuted form, according to the optional return
80 value @var{P}. For example, given the matrix @code{@var{A} = [1, 2; 3, 4]},
81 
82 @example
83 [@var{L}, @var{U}, @var{P}] = lu (@var{A})
84 @end example
85 
86 @noindent
87 returns
88 
89 @example
90 @group
91 L =
92 
93  1.00000 0.00000
94  0.33333 1.00000
95 
96 U =
97 
98  3.00000 4.00000
99  0.00000 0.66667
100 
101 P =
102 
103  0 1
104  1 0
105 @end group
106 @end example
107 
108 The matrix is not required to be square.
109 
110 When called with two or three output arguments and a sparse input matrix,
111 @code{lu} does not attempt to perform sparsity preserving column permutations.
112 Called with a fourth output argument, the sparsity preserving column
113 transformation @var{Q} is returned, such that
114 @code{@var{P} * @var{A} * @var{Q} = @var{L} * @var{U}}. This is the
115 @strong{preferred} way to call @code{lu} with sparse input matrices.
116 
117 Called with a fifth output argument and a sparse input matrix, @code{lu}
118 attempts to use a scaling factor @var{R} on the input matrix such that
119 @code{@var{P} * (@var{R} \ @var{A}) * @var{Q} = @var{L} * @var{U}}.
120 This typically leads to a sparser and more stable factorization.
121 
122 An additional input argument @var{thresh} that defines the pivoting
123 threshold can be given. @var{thresh} can be a scalar, in which case
124 it defines the @sc{umfpack} pivoting tolerance for both symmetric and
125 unsymmetric cases. If @var{thresh} is a 2-element vector, then the first
126 element defines the pivoting tolerance for the unsymmetric @sc{umfpack}
127 pivoting strategy and the second for the symmetric strategy. By default,
128 the values defined by @code{spparms} are used ([0.1, 0.001]).
129 
130 Given the string argument @qcode{"vector"}, @code{lu} returns the values
131 of @var{P} and @var{Q} as vector values, such that for full matrix,
132 @code{@var{A}(@var{P},:) = @var{L} * @var{U}}, and @code{@var{R}(@var{P},:)
133 * @var{A}(:,@var{Q}) = @var{L} * @var{U}}.
134 
135 With two output arguments, returns the permuted forms of the upper and
136 lower triangular matrices, such that @code{@var{A} = @var{L} * @var{U}}.
137 With one output argument @var{y}, then the matrix returned by the
138 @sc{lapack} routines is returned. If the input matrix is sparse then the
139 matrix @var{L} is embedded into @var{U} to give a return value similar to
140 the full case. For both full and sparse matrices, @code{lu} loses the
141 permutation information.
142 @seealso{luupdate, ilu, chol, hess, qr, qz, schur, svd}
143 @end deftypefn */)
144 {
145  int nargin = args.length ();
146  bool issparse = (nargin > 0 && args(0).issparse ());
147 
148  if (nargin < 1 || (issparse && nargin > 3) || (! issparse && nargin > 2))
149  print_usage ();
150 
151  bool vecout = false;
152  Matrix thresh;
153  int n = 1;
154 
155  while (n < nargin)
156  {
157  if (args(n).is_string ())
158  {
159  std::string tmp = args(n++).string_value ();
160 
161  if (tmp == "vector")
162  vecout = true;
163  else
164  error ("lu: unrecognized string argument");
165  }
166  else
167  {
168  if (! issparse)
169  error ("lu: can not define pivoting threshold THRESH for full matrices");
170 
171  Matrix tmp = args(n++).matrix_value ();
172  if (tmp.numel () == 1)
173  {
174  thresh.resize (1, 2);
175  thresh(0) = tmp(0);
176  thresh(1) = tmp(0);
177  }
178  else if (tmp.numel () == 2)
179  thresh = tmp;
180  else
181  error ("lu: THRESH must be a 1- or 2-element vector");
182  }
183  }
184 
185  octave_value_list retval;
186  bool scale = (nargout == 5);
187 
188  octave_value arg = args(0);
189 
190  octave_idx_type nr = arg.rows ();
191  octave_idx_type nc = arg.columns ();
192 
193  if (issparse)
194  {
195  if (arg.isempty ())
196  return octave_value_list (5, SparseMatrix ());
197 
198  if (arg.isreal ())
199  {
201 
202  if (nargout < 4)
203  {
204  warning_with_id ("Octave:lu:sparse_input",
205  "lu: function may fail when called with less than 4 output arguments and a sparse input");
206 
207  ColumnVector Qinit (nc);
208  for (octave_idx_type i = 0; i < nc; i++)
209  Qinit(i) = i;
210  math::sparse_lu<SparseMatrix> fact (m, Qinit, thresh, false,
211  true);
212 
213  if (nargout < 2)
214  retval(0) = fact.Y ();
215  else
216  {
217  retval.resize (nargout == 3 ? 3 : 2);
218  retval(1)
219  = octave_value (fact.U () * fact.Pc_mat ().transpose (),
221  nc, fact.col_perm ()));
222 
223  PermMatrix P = fact.Pr_mat ();
224  SparseMatrix L = fact.L ();
225 
226  if (nargout == 2)
227  retval(0)
228  = octave_value (P.transpose () * L,
230  nr, fact.row_perm ()));
231  else
232  {
233  retval(0) = L;
234  if (vecout)
235  retval(2) = fact.Pr_vec();
236  else
237  retval(2) = P;
238  }
239  }
240  }
241  else
242  {
243  retval.resize (scale ? 5 : 4);
244  math::sparse_lu<SparseMatrix> fact (m, thresh, scale);
245 
246  retval(0) = octave_value (fact.L (),
248  retval(1) = octave_value (fact.U (),
250 
251  if (vecout)
252  {
253  retval(2) = fact.Pr_vec ();
254  retval(3) = fact.Pc_vec ();
255  }
256  else
257  {
258  retval(2) = fact.Pr_mat ();
259  retval(3) = fact.Pc_mat ();
260  }
261 
262  if (scale)
263  retval(4) = fact.R ();
264  }
265  }
266  else if (arg.iscomplex ())
267  {
269 
270  if (nargout < 4)
271  {
272  warning_with_id ("Octave:lu:sparse_input",
273  "lu: function may fail when called with less than 4 output arguments and a sparse input");
274 
275  ColumnVector Qinit (nc);
276  for (octave_idx_type i = 0; i < nc; i++)
277  Qinit(i) = i;
278  math::sparse_lu<SparseComplexMatrix> fact (m, Qinit, thresh,
279  false, true);
280 
281  if (nargout < 2)
282  retval(0) = fact.Y ();
283  else
284  {
285  retval.resize (nargout == 3 ? 3 : 2);
286  retval(1)
287  = octave_value (fact.U () * fact.Pc_mat ().transpose (),
289  nc, fact.col_perm ()));
290 
291  PermMatrix P = fact.Pr_mat ();
292  SparseComplexMatrix L = fact.L ();
293  if (nargout == 2)
294  retval(0)
295  = octave_value (P.transpose () * L,
297  nr, fact.row_perm ()));
298  else
299  {
300  retval(0) = L;
301  if (vecout)
302  retval(2) = fact.Pr_vec();
303  else
304  retval(2) = P;
305  }
306  }
307  }
308  else
309  {
310  retval.resize (scale ? 5 : 4);
311  math::sparse_lu<SparseComplexMatrix> fact (m, thresh, scale);
312 
313  retval(0) = octave_value (fact.L (),
315  retval(1) = octave_value (fact.U (),
317 
318  if (vecout)
319  {
320  retval(2) = fact.Pr_vec ();
321  retval(3) = fact.Pc_vec ();
322  }
323  else
324  {
325  retval(2) = fact.Pr_mat ();
326  retval(3) = fact.Pc_mat ();
327  }
328 
329  if (scale)
330  retval(4) = fact.R ();
331  }
332 
333  }
334  else
335  err_wrong_type_arg ("lu", arg);
336  }
337  else
338  {
339  if (arg.isempty ())
340  return octave_value_list (3, Matrix ());
341 
342  if (arg.isreal ())
343  {
344  if (arg.is_single_type ())
345  {
347 
348  math::lu<FloatMatrix> fact (m);
349 
350  switch (nargout)
351  {
352  case 0:
353  case 1:
354  retval = ovl (fact.Y ());
355  break;
356 
357  case 2:
358  {
359  PermMatrix P = fact.P ();
360  FloatMatrix L = P.transpose () * fact.L ();
361  retval = ovl (L, get_lu_u (fact));
362  }
363  break;
364 
365  case 3:
366  default:
367  {
368  if (vecout)
369  retval = ovl (get_lu_l (fact), get_lu_u (fact),
370  fact.P_vec ());
371  else
372  retval = ovl (get_lu_l (fact), get_lu_u (fact),
373  fact.P ());
374  }
375  break;
376  }
377  }
378  else
379  {
380  Matrix m = arg.matrix_value ();
381 
382  math::lu<Matrix> fact (m);
383 
384  switch (nargout)
385  {
386  case 0:
387  case 1:
388  retval = ovl (fact.Y ());
389  break;
390 
391  case 2:
392  {
393  PermMatrix P = fact.P ();
394  Matrix L = P.transpose () * fact.L ();
395  retval = ovl (L, get_lu_u (fact));
396  }
397  break;
398 
399  case 3:
400  default:
401  {
402  if (vecout)
403  retval = ovl (get_lu_l (fact), get_lu_u (fact),
404  fact.P_vec ());
405  else
406  retval = ovl (get_lu_l (fact), get_lu_u (fact),
407  fact.P ());
408  }
409  break;
410  }
411  }
412  }
413  else if (arg.iscomplex ())
414  {
415  if (arg.is_single_type ())
416  {
418 
419  math::lu<FloatComplexMatrix> fact (m);
420 
421  switch (nargout)
422  {
423  case 0:
424  case 1:
425  retval = ovl (fact.Y ());
426  break;
427 
428  case 2:
429  {
430  PermMatrix P = fact.P ();
431  FloatComplexMatrix L = P.transpose () * fact.L ();
432  retval = ovl (L, get_lu_u (fact));
433  }
434  break;
435 
436  case 3:
437  default:
438  {
439  if (vecout)
440  retval = ovl (get_lu_l (fact), get_lu_u (fact),
441  fact.P_vec ());
442  else
443  retval = ovl (get_lu_l (fact), get_lu_u (fact),
444  fact.P ());
445  }
446  break;
447  }
448  }
449  else
450  {
452 
453  math::lu<ComplexMatrix> fact (m);
454 
455  switch (nargout)
456  {
457  case 0:
458  case 1:
459  retval = ovl (fact.Y ());
460  break;
461 
462  case 2:
463  {
464  PermMatrix P = fact.P ();
465  ComplexMatrix L = P.transpose () * fact.L ();
466  retval = ovl (L, get_lu_u (fact));
467  }
468  break;
469 
470  case 3:
471  default:
472  {
473  if (vecout)
474  retval = ovl (get_lu_l (fact), get_lu_u (fact),
475  fact.P_vec ());
476  else
477  retval = ovl (get_lu_l (fact), get_lu_u (fact),
478  fact.P ());
479  }
480  break;
481  }
482  }
483  }
484  else
485  err_wrong_type_arg ("lu", arg);
486  }
487 
488  return retval;
489 }
490 
491 /*
492 %!assert(lu ([1, 2; 3, 4]), [3, 4; 1/3, 2/3], eps)
493 
494 %!test
495 %! [l, u] = lu ([1, 2; 3, 4]);
496 %! assert (l, [1/3, 1; 1, 0], sqrt (eps));
497 %! assert (u, [3, 4; 0, 2/3], sqrt (eps));
498 
499 %!test
500 %! [l, u, p] = lu ([1, 2; 3, 4]);
501 %! assert (l, [1, 0; 1/3, 1], sqrt (eps));
502 %! assert (u, [3, 4; 0, 2/3], sqrt (eps));
503 %! assert (p(:,:), [0, 1; 1, 0], sqrt (eps));
504 
505 %!test
506 %! [l, u, p] = lu ([1, 2; 3, 4], "vector");
507 %! assert (l, [1, 0; 1/3, 1], sqrt (eps));
508 %! assert (u, [3, 4; 0, 2/3], sqrt (eps));
509 %! assert (p, [2;1], sqrt (eps));
510 
511 %!test
512 %! [l, u, p] = lu ([1, 2; 3, 4; 5, 6]);
513 %! assert (l, [1, 0; 1/5, 1; 3/5, 1/2], sqrt (eps));
514 %! assert (u, [5, 6; 0, 4/5], sqrt (eps));
515 %! assert (p(:,:), [0, 0, 1; 1, 0, 0; 0 1 0], sqrt (eps));
516 
517 %!assert (lu (single ([1, 2; 3, 4])), single ([3, 4; 1/3, 2/3]), eps ("single"))
518 
519 %!test
520 %! [l, u] = lu (single ([1, 2; 3, 4]));
521 %! assert (l, single ([1/3, 1; 1, 0]), sqrt (eps ("single")));
522 %! assert (u, single ([3, 4; 0, 2/3]), sqrt (eps ("single")));
523 
524 %!test
525 %! [l, u, p] = lu (single ([1, 2; 3, 4]));
526 %! assert (l, single ([1, 0; 1/3, 1]), sqrt (eps ("single")));
527 %! assert (u, single ([3, 4; 0, 2/3]), sqrt (eps ("single")));
528 %! assert (p(:,:), single ([0, 1; 1, 0]), sqrt (eps ("single")));
529 
530 %!test
531 %! [l, u, p] = lu (single ([1, 2; 3, 4]), "vector");
532 %! assert (l, single ([1, 0; 1/3, 1]), sqrt (eps ("single")));
533 %! assert (u, single ([3, 4; 0, 2/3]), sqrt (eps ("single")));
534 %! assert (p, single ([2;1]), sqrt (eps ("single")));
535 
536 %!test
537 %! [l u p] = lu (single ([1, 2; 3, 4; 5, 6]));
538 %! assert (l, single ([1, 0; 1/5, 1; 3/5, 1/2]), sqrt (eps ("single")));
539 %! assert (u, single ([5, 6; 0, 4/5]), sqrt (eps ("single")));
540 %! assert (p(:,:), single ([0, 0, 1; 1, 0, 0; 0 1 0]), sqrt (eps ("single")));
541 
542 # complex matrix input
543 %!test
544 %! C = [1, 0, 1, 0;
545 %! 1i, 1/3, -1i, 1/3;
546 %! 1, 2i/3, 1, -2i/3;
547 %! 1i, -1/3, -1i, -1/3];
548 %! [L, U, P] = lu (C);
549 %! assert (rcond (C), 1/8, eps);
550 %! assert (norm (P'*L*U - C, Inf) < eps);
551 
552 %!testif HAVE_UMFPACK
553 %! Bi = [1 2 3 4 5 2 3 6 7 8 4 5 7 8 9];
554 %! Bj = [1 3 4 5 6 7 8 9 11 12 13 14 15 16 17];
555 %! Bv = [1 1 1 1 1 1 -1 1 1 1 1 -1 1 -1 1];
556 %! B = sparse (Bi, Bj, Bv);
557 %! warning ("off", "Octave:lu:sparse_input", "local");
558 %! [L, U] = lu (B);
559 %! assert (L*U, B);
560 %! [L, U, P] = lu(B);
561 %! assert (P'*L*U, B);
562 %! [L, U, P, Q] = lu (B);
563 %! assert (P'*L*U*Q', B);
564 
565 %!error lu ()
566 %!testif HAVE_UMFPACK
567 %! fail ("[l,u] = lu (sparse (magic (3)))", "warning", "function may fail");
568 %!error <can not define pivoting threshold> lu ([1, 2; 3, 4], 2)
569 
570 */
571 
572 static bool
573 check_lu_dims (const octave_value& l, const octave_value& u,
574  const octave_value& p)
575 {
576  octave_idx_type m = l.rows ();
577  octave_idx_type k = u.rows ();
578  octave_idx_type n = u.columns ();
579 
580  return ((l.ndims () == 2 && u.ndims () == 2 && k == l.columns ())
581  && k == std::min (m, n)
582  && (p.is_undefined () || p.rows () == m));
583 }
584 
585 DEFUN (luupdate, args, ,
586  doc: /* -*- texinfo -*-
587 @deftypefn {} {[@var{L}, @var{U}] =} luupdate (@var{L}, @var{U}, @var{x}, @var{y})
588 @deftypefnx {} {[@var{L}, @var{U}, @var{P}] =} luupdate (@var{L}, @var{U}, @var{P}, @var{x}, @var{y})
589 Given an LU@tie{}factorization of a real or complex matrix
590 @w{@var{A} = @var{L}*@var{U}}, @var{L}@tie{}lower unit trapezoidal and
591 @var{U}@tie{}upper trapezoidal, return the LU@tie{}factorization
592 of @w{@var{A} + @var{x}*@var{y}.'}, where @var{x} and @var{y} are
593 column vectors (rank-1 update) or matrices with equal number of columns
594 (rank-k update).
595 
596 Optionally, row-pivoted updating can be used by supplying a row permutation
597 (pivoting) matrix @var{P}; in that case, an updated permutation matrix is
598 returned. Note that if @var{L}, @var{U}, @var{P} is a pivoted
599 LU@tie{}factorization as obtained by @code{lu}:
600 
601 @example
602 [@var{L}, @var{U}, @var{P}] = lu (@var{A});
603 @end example
604 
605 @noindent
606 then a factorization of @tcode{@var{A}+@var{x}*@var{y}.'} can be obtained
607 either as
608 
609 @example
610 [@var{L1}, @var{U1}] = lu (@var{L}, @var{U}, @var{P}*@var{x}, @var{y})
611 @end example
612 
613 @noindent
614 or
615 
616 @example
617 [@var{L1}, @var{U1}, @var{P1}] = lu (@var{L}, @var{U}, @var{P}, @var{x}, @var{y})
618 @end example
619 
620 The first form uses the unpivoted algorithm, which is faster, but less
621 stable. The second form uses a slower pivoted algorithm, which is more
622 stable.
623 
624 The matrix case is done as a sequence of rank-1 updates; thus, for large
625 enough k, it will be both faster and more accurate to recompute the
626 factorization from scratch.
627 @seealso{lu, cholupdate, qrupdate}
628 @end deftypefn */)
629 {
630  int nargin = args.length ();
631 
632  if (nargin != 4 && nargin != 5)
633  print_usage ();
634 
635  bool pivoted = (nargin == 5);
636 
637  octave_value argl = args(0);
638  octave_value argu = args(1);
639  octave_value argp = (pivoted ? args(2) : octave_value ());
640  octave_value argx = args(2 + pivoted);
641  octave_value argy = args(3 + pivoted);
642 
643  if (! (argl.isnumeric () && argu.isnumeric ()
644  && argx.isnumeric () && argy.isnumeric ()
645  && (! pivoted || argp.is_perm_matrix ())))
646  error ("luupdate: L, U, X, and Y must be numeric");
647 
648  if (! check_lu_dims (argl, argu, argp))
649  error ("luupdate: dimension mismatch");
650 
651  PermMatrix P = (pivoted
652  ? argp.perm_matrix_value ()
653  : PermMatrix::eye (argl.rows ()));
654 
655  if (argl.isreal () && argu.isreal ()
656  && argx.isreal () && argy.isreal ())
657  {
658  // all real case
659  if (argl.is_single_type () || argu.is_single_type ()
660  || argx.is_single_type () || argy.is_single_type ())
661  {
662  FloatMatrix L = argl.float_matrix_value ();
663  FloatMatrix U = argu.float_matrix_value ();
664  FloatMatrix x = argx.float_matrix_value ();
665  FloatMatrix y = argy.float_matrix_value ();
666 
667  math::lu<FloatMatrix> fact (L, U, P);
668  if (pivoted)
669  fact.update_piv (x, y);
670  else
671  fact.update (x, y);
672 
673  if (pivoted)
674  return ovl (get_lu_l (fact), get_lu_u (fact), fact.P ());
675  else
676  return ovl (get_lu_l (fact), get_lu_u (fact));
677  }
678  else
679  {
680  Matrix L = argl.matrix_value ();
681  Matrix U = argu.matrix_value ();
682  Matrix x = argx.matrix_value ();
683  Matrix y = argy.matrix_value ();
684 
685  math::lu<Matrix> fact (L, U, P);
686  if (pivoted)
687  fact.update_piv (x, y);
688  else
689  fact.update (x, y);
690 
691  if (pivoted)
692  return ovl (get_lu_l (fact), get_lu_u (fact), fact.P ());
693  else
694  return ovl (get_lu_l (fact), get_lu_u (fact));
695  }
696  }
697  else
698  {
699  // complex case
700  if (argl.is_single_type () || argu.is_single_type ()
701  || argx.is_single_type () || argy.is_single_type ())
702  {
707 
708  math::lu<FloatComplexMatrix> fact (L, U, P);
709  if (pivoted)
710  fact.update_piv (x, y);
711  else
712  fact.update (x, y);
713 
714  if (pivoted)
715  return ovl (get_lu_l (fact), get_lu_u (fact), fact.P ());
716  else
717  return ovl (get_lu_l (fact), get_lu_u (fact));
718  }
719  else
720  {
725 
726  math::lu<ComplexMatrix> fact (L, U, P);
727  if (pivoted)
728  fact.update_piv (x, y);
729  else
730  fact.update (x, y);
731 
732  if (pivoted)
733  return ovl (get_lu_l (fact), get_lu_u (fact), fact.P ());
734  else
735  return ovl (get_lu_l (fact), get_lu_u (fact));
736  }
737  }
738 }
739 
740 /*
741 %!shared A, u, v, Ac, uc, vc
742 %! A = [0.091364 0.613038 0.999083;
743 %! 0.594638 0.425302 0.603537;
744 %! 0.383594 0.291238 0.085574;
745 %! 0.265712 0.268003 0.238409;
746 %! 0.669966 0.743851 0.445057 ];
747 %!
748 %! u = [0.85082;
749 %! 0.76426;
750 %! 0.42883;
751 %! 0.53010;
752 %! 0.80683 ];
753 %!
754 %! v = [0.98810;
755 %! 0.24295;
756 %! 0.43167 ];
757 %!
758 %! Ac = [0.620405 + 0.956953i 0.480013 + 0.048806i 0.402627 + 0.338171i;
759 %! 0.589077 + 0.658457i 0.013205 + 0.279323i 0.229284 + 0.721929i;
760 %! 0.092758 + 0.345687i 0.928679 + 0.241052i 0.764536 + 0.832406i;
761 %! 0.912098 + 0.721024i 0.049018 + 0.269452i 0.730029 + 0.796517i;
762 %! 0.112849 + 0.603871i 0.486352 + 0.142337i 0.355646 + 0.151496i ];
763 %!
764 %! uc = [0.20351 + 0.05401i;
765 %! 0.13141 + 0.43708i;
766 %! 0.29808 + 0.08789i;
767 %! 0.69821 + 0.38844i;
768 %! 0.74871 + 0.25821i ];
769 %!
770 %! vc = [0.85839 + 0.29468i;
771 %! 0.20820 + 0.93090i;
772 %! 0.86184 + 0.34689i ];
773 %!
774 
775 %!testif HAVE_QRUPDATE_LUU
776 %! [L,U,P] = lu (A);
777 %! [L,U] = luupdate (L,U,P*u,v);
778 %! assert (norm (vec (tril (L)-L), Inf) == 0);
779 %! assert (norm (vec (triu (U)-U), Inf) == 0);
780 %! assert (norm (vec (P'*L*U - A - u*v.'), Inf) < norm (A)*1e1*eps);
781 %!
782 %!testif HAVE_QRUPDATE_LUU
783 %! [L,U,P] = lu (Ac);
784 %! [L,U] = luupdate (L,U,P*uc,vc);
785 %! assert (norm (vec (tril (L)-L), Inf) == 0);
786 %! assert (norm (vec (triu (U)-U), Inf) == 0);
787 %! assert (norm (vec (P'*L*U - Ac - uc*vc.'), Inf) < norm (Ac)*1e1*eps);
788 
789 %!testif HAVE_QRUPDATE_LUU
790 %! [L,U,P] = lu (single (A));
791 %! [L,U] = luupdate (L,U,P* single (u), single (v));
792 %! assert (norm (vec (tril (L)-L), Inf) == 0);
793 %! assert (norm (vec (triu (U)-U), Inf) == 0);
794 %! assert (norm (vec (P'*L*U - single (A) - single (u)* single (v).'), Inf)
795 %! < norm (single (A))*1e1* eps ("single"));
796 %!
797 %!testif HAVE_QRUPDATE_LUU
798 %! [L,U,P] = lu (single (Ac));
799 %! [L,U] = luupdate (L,U,P* single (uc),single (vc));
800 %! assert (norm (vec (tril (L)-L), Inf) == 0);
801 %! assert (norm (vec (triu (U)-U), Inf) == 0);
802 %! assert (norm (vec (P'*L*U - single (Ac) - single (uc)* single (vc).'), Inf)
803 %! < norm (single (Ac))*1e1* eps ("single"));
804 
805 %!testif HAVE_QRUPDATE_LUU
806 %! [L,U,P] = lu (A);
807 %! [L,U,P] = luupdate (L,U,P,u,v);
808 %! assert (norm (vec (tril (L)-L), Inf) == 0);
809 %! assert (norm (vec (triu (U)-U), Inf) == 0);
810 %! assert (norm (vec (P'*L*U - A - u*v.'), Inf) < norm (A)*1e1*eps);
811 %!
812 %!testif HAVE_QRUPDATE_LUU
813 %! [L,U,P] = lu (A);
814 %! [~,ordcols] = max (P,[],1);
815 %! [~,ordrows] = max (P,[],2);
816 %! P1 = eye (size (P))(:,ordcols);
817 %! P2 = eye (size (P))(ordrows,:);
818 %! assert (P1 == P);
819 %! assert (P2 == P);
820 %! [L,U,P] = luupdate (L,U,P,u,v);
821 %! [L,U,P1] = luupdate (L,U,P1,u,v);
822 %! [L,U,P2] = luupdate (L,U,P2,u,v);
823 %! assert (P1 == P);
824 %! assert (P2 == P);
825 %!
826 %!testif HAVE_QRUPDATE_LUU
827 %! [L,U,P] = lu (Ac);
828 %! [L,U,P] = luupdate (L,U,P,uc,vc);
829 %! assert (norm (vec (tril (L)-L), Inf) == 0);
830 %! assert (norm (vec (triu (U)-U), Inf) == 0);
831 %! assert (norm (vec (P'*L*U - Ac - uc*vc.'), Inf) < norm (Ac)*1e1*eps);
832 
833 %!testif HAVE_QRUPDATE_LUU
834 %! [L,U,P] = lu (single (A));
835 %! [L,U,P] = luupdate (L,U,P,single (u),single (v));
836 %! assert (norm (vec (tril (L)-L), Inf) == 0);
837 %! assert (norm (vec (triu (U)-U), Inf) == 0);
838 %! assert (norm (vec (P'*L*U - single (A) - single (u)* single (v).'), Inf)
839 %! < norm (single (A))*1e1* eps ("single"));
840 %!
841 %!testif HAVE_QRUPDATE_LUU
842 %! [L,U,P] = lu (single (Ac));
843 %! [L,U,P] = luupdate (L,U,P,single (uc),single (vc));
844 %! assert (norm (vec (tril (L)-L), Inf) == 0);
845 %! assert (norm (vec (triu (U)-U), Inf) == 0);
846 %! assert (norm (vec (P'*L*U - single (Ac) - single (uc)* single (vc).'), Inf)
847 %! < norm (single (Ac))*1e1* eps ("single"));
848 */
849 
850 OCTAVE_END_NAMESPACE(octave)
charNDArray min(char d, const charNDArray &m)
Definition: chNDArray.cc:207
octave_idx_type numel() const
Number of elements in the array.
Definition: Array.h:414
@ Permuted_Lower
Definition: MatrixType.h:48
@ Permuted_Upper
Definition: MatrixType.h:47
Definition: dMatrix.h:42
void resize(octave_idx_type nr, octave_idx_type nc, double rfv=0)
Definition: dMatrix.h:158
static PermMatrix eye(octave_idx_type n)
Definition: PermMatrix.cc:205
PermMatrix transpose() const
Definition: PermMatrix.cc:101
Definition: lu.h:42
void resize(octave_idx_type n, const octave_value &rfv=octave_value())
Definition: ovl.h:117
bool is_undefined() const
Definition: ov.h:595
SparseMatrix sparse_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:900
int ndims() const
Definition: ov.h:551
octave_idx_type rows() const
Definition: ov.h:545
bool isreal() const
Definition: ov.h:738
bool isnumeric() const
Definition: ov.h:750
ComplexMatrix complex_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:871
bool is_perm_matrix() const
Definition: ov.h:634
bool is_single_type() const
Definition: ov.h:698
bool isempty() const
Definition: ov.h:601
FloatMatrix float_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:856
bool iscomplex() const
Definition: ov.h:741
octave_idx_type columns() const
Definition: ov.h:547
Matrix matrix_value(bool frc_str_conv=false) const
Definition: ov.h:853
FloatComplexMatrix float_complex_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:875
PermMatrix perm_matrix_value() const
Definition: ov.h:923
SparseComplexMatrix sparse_complex_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:904
OCTAVE_BEGIN_NAMESPACE(octave) static octave_value daspk_fcn
void print_usage(void)
Definition: defun-int.h:72
#define DEFUN(name, args_name, nargout_name, doc)
Macro to define a builtin function.
Definition: defun.h:56
void warning_with_id(const char *id, const char *fmt,...)
Definition: error.cc:1078
void() error(const char *fmt,...)
Definition: error.cc:988
void err_wrong_type_arg(const char *name, const char *s)
Definition: errwarn.cc:166
void scale(Matrix &m, double x, double y, double z)
Definition: graphics.cc:5500
F77_RET_T const F77_DBLE * x
T octave_idx_type m
Definition: mx-inlines.cc:781
octave_idx_type n
Definition: mx-inlines.cc:761
return octave_value(v1.char_array_value() . concat(v2.char_array_value(), ra_idx),((a1.is_sq_string()||a2.is_sq_string()) ? '\'' :'"'))
octave_value_list ovl(const OV_Args &... args)
Construct an octave_value_list with less typing.
Definition: ovl.h:219