GNU Octave  9.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
chol.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 <string>
31 
32 #include "Matrix.h"
33 #include "chol.h"
34 #include "oct-string.h"
35 #include "sparse-chol.h"
36 #include "sparse-util.h"
37 
38 #include "defun.h"
39 #include "error.h"
40 #include "errwarn.h"
41 #include "ov.h"
42 #include "ovl.h"
43 
45 
46 template <typename CHOLT>
47 static octave_value
48 get_chol (const CHOLT& fact)
49 {
50  return octave_value (fact.chol_matrix());
51 }
52 
53 template <typename CHOLT>
54 static octave_value
55 get_chol_r (const CHOLT& fact)
56 {
57  return octave_value (fact.chol_matrix (),
59 }
60 
61 template <typename CHOLT>
62 static octave_value
63 get_chol_l (const CHOLT& fact)
64 {
65  return octave_value (fact.chol_matrix ().transpose (),
67 }
68 
69 DEFUN (chol, args, nargout,
70  doc: /* -*- texinfo -*-
71 @deftypefn {} {@var{R} =} chol (@var{A})
72 @deftypefnx {} {[@var{R}, @var{p}] =} chol (@var{A})
73 @deftypefnx {} {[@var{R}, @var{p}, @var{Q}] =} chol (@var{A})
74 @deftypefnx {} {[@var{R}, @var{p}, @var{Q}] =} chol (@var{A}, "vector")
75 @deftypefnx {} {[@var{L}, @dots{}] =} chol (@dots{}, "lower")
76 @deftypefnx {} {[@var{R}, @dots{}] =} chol (@dots{}, "upper")
77 @cindex Cholesky factorization
78 Compute the upper Cholesky@tie{}factor, @var{R}, of the real symmetric
79 or complex Hermitian positive definite matrix @var{A}.
80 
81 The upper Cholesky@tie{}factor @var{R} is computed by using the upper
82 triangular part of matrix @var{A} and is defined by
83 @tex
84 $ R^T R = A $.
85 @end tex
86 @ifnottex
87 
88 @example
89 @var{R}' * @var{R} = @var{A}.
90 @end example
91 
92 @end ifnottex
93 
94 Calling @code{chol} using the optional @qcode{"upper"} flag has the
95 same behavior. In contrast, using the optional @qcode{"lower"} flag,
96 @code{chol} returns the lower triangular factorization, computed by using
97 the lower triangular part of matrix @var{A}, such that
98 @tex
99 $ L L^T = A $.
100 @end tex
101 @ifnottex
102 
103 @example
104 @var{L} * @var{L}' = @var{A}.
105 @end example
106 
107 @end ifnottex
108 
109 Called with one output argument @code{chol} fails if matrix @var{A} is
110 not positive definite. Note that if matrix @var{A} is not real symmetric
111 or complex Hermitian then the lower triangular part is considered to be
112 the (complex conjugate) transpose of the upper triangular part, or vice
113 versa, given the @qcode{"lower"} flag.
114 
115 Called with two or more output arguments @var{p} flags whether the matrix
116 @var{A} was positive definite and @code{chol} does not fail. A zero value
117 of @var{p} indicates that matrix @var{A} is positive definite and @var{R}
118 gives the factorization. Otherwise, @var{p} will have a positive value.
119 
120 If called with three output arguments matrix @var{A} must be sparse and
121 a sparsity preserving row/column permutation is applied to matrix @var{A}
122 prior to the factorization. That is @var{R} is the factorization of
123 @code{@var{A}(@var{Q},@var{Q})} such that
124 @tex
125 $ R^T R = Q^T A Q$.
126 @end tex
127 @ifnottex
128 
129 @example
130 @var{R}' * @var{R} = @var{Q}' * @var{A} * @var{Q}.
131 @end example
132 
133 @end ifnottex
134 
135 The sparsity preserving permutation is generally returned as a matrix.
136 However, given the optional flag @qcode{"vector"}, @var{Q} will be
137 returned as a vector such that
138 @tex
139 $ R^T R = A (Q, Q)$.
140 @end tex
141 @ifnottex
142 
143 @example
144 @var{R}' * @var{R} = @var{A}(@var{Q}, @var{Q}).
145 @end example
146 
147 @end ifnottex
148 
149 In general the lower triangular factorization is significantly faster for
150 sparse matrices.
151 @seealso{hess, lu, qr, qz, schur, svd, ichol, cholinv, chol2inv, cholupdate,
152 cholinsert, choldelete, cholshift}
153 @end deftypefn */)
154 {
155  int nargin = args.length ();
156 
157  if (nargin < 1 || nargin > 3 || nargout > 3)
158  print_usage ();
159  if (nargout > 2 && ! args(0).issparse ())
160  error ("chol: using three output arguments, matrix A must be sparse");
161 
162  bool LLt = false;
163  bool vecout = false;
164 
165  int n = 1;
166  while (n < nargin)
167  {
168  std::string tmp = args(n++).xstring_value ("chol: optional arguments must be strings");
169 
170  if (string::strcmpi (tmp, "vector"))
171  vecout = true;
172  else if (string::strcmpi (tmp, "lower"))
173  LLt = true;
174  else if (string::strcmpi (tmp, "upper"))
175  LLt = false;
176  else
177  error (R"(chol: optional argument must be one of "vector", "lower", or "upper")");
178  }
179 
180  octave_value_list retval;
181  octave_value arg = args(0);
182 
183  if (arg.isempty ())
184  return ovl (Matrix ());
185 
186  if (arg.issparse ())
187  {
188  octave_idx_type info;
189  bool natural = (nargout != 3);
190  bool force = nargout > 1;
191 
192  if (arg.isreal ())
193  {
195 
196  math::sparse_chol<SparseMatrix> fact (m, info, natural, force);
197 
198  if (nargout == 3)
199  {
200  if (vecout)
201  retval(2) = fact.perm ();
202  else
203  retval(2) = fact.Q ();
204  }
205 
206  if (nargout >= 2 || info == 0)
207  {
208  retval(1) = info;
209  if (LLt)
210  retval(0) = fact.L ();
211  else
212  retval(0) = fact.R ();
213  }
214  else
215  error ("chol: input matrix must be positive definite");
216  }
217  else if (arg.iscomplex ())
218  {
220 
221  math::sparse_chol<SparseComplexMatrix> fact (m, info, natural, force);
222 
223  if (nargout == 3)
224  {
225  if (vecout)
226  retval(2) = fact.perm ();
227  else
228  retval(2) = fact.Q ();
229  }
230 
231  if (nargout >= 2 || info == 0)
232  {
233  retval(1) = info;
234  if (LLt)
235  retval(0) = fact.L ();
236  else
237  retval(0) = fact.R ();
238  }
239  else
240  error ("chol: input matrix must be positive definite");
241  }
242  else
243  err_wrong_type_arg ("chol", arg);
244  }
245  else if (arg.is_single_type ())
246  {
247  if (vecout)
248  error (R"(chol: A must be sparse for the "vector" option)");
249  if (arg.isreal ())
250  {
252 
253  octave_idx_type info;
254 
255  math::chol<FloatMatrix> fact (m, info, LLt != true);
256 
257  if (nargout == 2 || info == 0)
258  retval = ovl (get_chol (fact), info);
259  else
260  error ("chol: input matrix must be positive definite");
261  }
262  else if (arg.iscomplex ())
263  {
265 
266  octave_idx_type info;
267 
268  math::chol<FloatComplexMatrix> fact (m, info, LLt != true);
269 
270  if (nargout == 2 || info == 0)
271  retval = ovl (get_chol (fact), info);
272  else
273  error ("chol: input matrix must be positive definite");
274  }
275  else
276  err_wrong_type_arg ("chol", arg);
277  }
278  else
279  {
280  if (vecout)
281  error (R"(chol: A must be sparse for the "vector" option)");
282  if (arg.isreal ())
283  {
284  Matrix m = arg.matrix_value ();
285 
286  octave_idx_type info;
287 
288  math::chol<Matrix> fact (m, info, LLt != true);
289 
290  if (nargout == 2 || info == 0)
291  retval = ovl (get_chol (fact), info);
292  else
293  error ("chol: input matrix must be positive definite");
294  }
295  else if (arg.iscomplex ())
296  {
298 
299  octave_idx_type info;
300 
301  math::chol<ComplexMatrix> fact (m, info, LLt != true);
302 
303  if (nargout == 2 || info == 0)
304  retval = ovl (get_chol (fact), info);
305  else
306  error ("chol: input matrix must be positive definite");
307  }
308  else
309  err_wrong_type_arg ("chol", arg);
310  }
311 
312  return retval;
313 }
314 
315 /*
316 %!assert (chol ([2, 1; 1, 1]), [sqrt(2), 1/sqrt(2); 0, 1/sqrt(2)], sqrt (eps))
317 %!assert (chol (single ([2, 1; 1, 1])),
318 %! single ([sqrt(2), 1/sqrt(2); 0, 1/sqrt(2)]), sqrt (eps ("single")))
319 
320 %!assert (chol ([2, 1; 1, 1], "upper"), [sqrt(2), 1/sqrt(2); 0, 1/sqrt(2)],
321 %! sqrt (eps))
322 %!assert (chol ([2, 1; 1, 1], "lower"), [sqrt(2), 0; 1/sqrt(2), 1/sqrt(2)],
323 %! sqrt (eps))
324 
325 %!assert (chol ([2, 1; 1, 1], "lower"), chol ([2, 1; 1, 1], "LoweR"))
326 %!assert (chol ([2, 1; 1, 1], "upper"), chol ([2, 1; 1, 1], "Upper"))
327 
328 ## Check the "vector" option which only affects the 3rd argument and
329 ## is only valid for sparse input.
330 %!testif HAVE_CHOLMOD
331 %! a = sparse ([2 1; 1 1]);
332 %! r = sparse ([sqrt(2), 1/sqrt(2); 0, 1/sqrt(2)]);
333 %! [rd, pd, qd] = chol (a);
334 %! [rv, pv, qv] = chol (a, "vector");
335 %! assert (r, rd, eps)
336 %! assert (r, rv, eps)
337 %! assert (pd, 0)
338 %! assert (pd, pv)
339 %! assert (qd, sparse (eye (2)))
340 %! assert (qv, [1 2])
341 %!
342 %! [rv, pv, qv] = chol (a, "Vector"); # check case sensitivity
343 %! assert (r, rv, eps)
344 %! assert (pd, pv)
345 %! assert (qv, [1 2])
346 
347 %!testif HAVE_CHOLMOD <*42587>
348 %! A = sparse ([1 0 8;0 1 8;8 8 1]);
349 %! [Q, p] = chol (A);
350 %! assert (p != 0);
351 
352 %!error chol ()
353 %!error <matrix must be positive definite> chol ([1, 2; 3, 4])
354 %!error <requires square matrix> chol ([1, 2; 3, 4; 5, 6])
355 %!error <optional arguments must be strings> chol (1, 2)
356 %!error <optional argument must be one of "vector", "lower"> chol (1, "foobar")
357 %!error <matrix A must be sparse> [L,p,Q] = chol ([1, 2; 3, 4])
358 %!error <A must be sparse> [L, p] = chol ([1, 2; 3, 4], "vector")
359 */
360 
361 DEFUN (cholinv, args, ,
362  doc: /* -*- texinfo -*-
363 @deftypefn {} {@var{Ainv} =} cholinv (@var{A})
364 Compute the inverse of the symmetric positive definite matrix @var{A} using
365 the Cholesky@tie{}factorization.
366 @seealso{chol, chol2inv, inv}
367 @end deftypefn */)
368 {
369  if (args.length () != 1)
370  print_usage ();
371 
372  octave_value retval;
373  octave_value arg = args(0);
374 
375  octave_idx_type nr = arg.rows ();
376  octave_idx_type nc = arg.columns ();
377 
378  if (nr == 0 || nc == 0)
379  retval = Matrix ();
380  else
381  {
382  if (arg.issparse ())
383  {
384  octave_idx_type info;
385 
386  if (arg.isreal ())
387  {
389 
390  math::sparse_chol<SparseMatrix> chol (m, info);
391 
392  if (info == 0)
393  retval = chol.inverse ();
394  else
395  error ("cholinv: A must be positive definite");
396  }
397  else if (arg.iscomplex ())
398  {
400 
401  math::sparse_chol<SparseComplexMatrix> chol (m, info);
402 
403  if (info == 0)
404  retval = chol.inverse ();
405  else
406  error ("cholinv: A must be positive definite");
407  }
408  else
409  err_wrong_type_arg ("cholinv", arg);
410  }
411  else if (arg.is_single_type ())
412  {
413  if (arg.isreal ())
414  {
416 
417  octave_idx_type info;
418  math::chol<FloatMatrix> chol (m, info);
419  if (info == 0)
420  retval = chol.inverse ();
421  else
422  error ("cholinv: A must be positive definite");
423  }
424  else if (arg.iscomplex ())
425  {
427 
428  octave_idx_type info;
429  math::chol<FloatComplexMatrix> chol (m, info);
430  if (info == 0)
431  retval = chol.inverse ();
432  else
433  error ("cholinv: A must be positive definite");
434  }
435  else
436  err_wrong_type_arg ("chol", arg);
437  }
438  else
439  {
440  if (arg.isreal ())
441  {
442  Matrix m = arg.matrix_value ();
443 
444  octave_idx_type info;
445  math::chol<Matrix> chol (m, info);
446  if (info == 0)
447  retval = chol.inverse ();
448  else
449  error ("cholinv: A must be positive definite");
450  }
451  else if (arg.iscomplex ())
452  {
454 
455  octave_idx_type info;
456  math::chol<ComplexMatrix> chol (m, info);
457  if (info == 0)
458  retval = chol.inverse ();
459  else
460  error ("cholinv: A must be positive definite");
461  }
462  else
463  err_wrong_type_arg ("chol", arg);
464  }
465  }
466 
467  return retval;
468 }
469 
470 /*
471 %!shared A, Ainv
472 %! A = [2,0.2;0.2,1];
473 %! Ainv = inv (A);
474 %!test
475 %! Ainv1 = cholinv (A);
476 %! assert (norm (Ainv-Ainv1), 0, 1e-10);
477 %!testif HAVE_CHOLMOD
478 %! Ainv2 = inv (sparse (A));
479 %! assert (norm (Ainv-Ainv2), 0, 1e-10);
480 %!testif HAVE_CHOLMOD
481 %! Ainv3 = cholinv (sparse (A));
482 %! assert (norm (Ainv-Ainv3), 0, 1e-10);
483 */
484 
485 DEFUN (chol2inv, args, ,
486  doc: /* -*- texinfo -*-
487 @deftypefn {} {@var{Ainv} =} chol2inv (@var{R})
488 Invert a symmetric, positive definite square matrix from its Cholesky
489 decomposition, @var{R}.
490 
491 Note that @var{R} should be an upper-triangular matrix with positive diagonal
492 elements. @code{chol2inv (@var{U})} provides @code{inv (@var{R}'*@var{R})} but
493 is much faster than using @code{inv}.
494 @seealso{chol, cholinv, inv}
495 @end deftypefn */)
496 {
497  if (args.length () != 1)
498  print_usage ();
499 
500  octave_value retval;
501 
502  octave_value arg = args(0);
503 
504  octave_idx_type nr = arg.rows ();
505  octave_idx_type nc = arg.columns ();
506 
507  if (nr == 0 || nc == 0)
508  retval = Matrix ();
509  else
510  {
511  if (arg.issparse ())
512  {
513  if (arg.isreal ())
514  {
516 
517  retval = math::chol2inv (r);
518  }
519  else if (arg.iscomplex ())
520  {
522 
523  retval = math::chol2inv (r);
524  }
525  else
526  err_wrong_type_arg ("chol2inv", arg);
527  }
528  else if (arg.is_single_type ())
529  {
530  if (arg.isreal ())
531  {
533 
534  retval = math::chol2inv (r);
535  }
536  else if (arg.iscomplex ())
537  {
539 
540  retval = math::chol2inv (r);
541  }
542  else
543  err_wrong_type_arg ("chol2inv", arg);
544 
545  }
546  else
547  {
548  if (arg.isreal ())
549  {
550  Matrix r = arg.matrix_value ();
551 
552  retval = math::chol2inv (r);
553  }
554  else if (arg.iscomplex ())
555  {
557 
558  retval = math::chol2inv (r);
559  }
560  else
561  err_wrong_type_arg ("chol2inv", arg);
562  }
563  }
564 
565  return retval;
566 }
567 
568 /*
569 
570 ## Test for bug #36437
571 %!function sparse_chol2inv (T, tol)
572 %! iT = inv (T);
573 %! ciT = chol2inv (chol (T));
574 %! assert (ciT, iT, tol);
575 %! assert (chol2inv (chol ( full (T))), ciT, tol*2);
576 %!endfunction
577 
578 %!testif HAVE_CHOLMOD
579 %! A = gallery ("poisson", 3);
580 %! sparse_chol2inv (A, eps);
581 
582 %!testif HAVE_CHOLMOD
583 %! n = 10;
584 %! B = spdiags (ones (n, 1) * [1 2 1], [-1 0 1], n, n);
585 %! sparse_chol2inv (B, eps*100);
586 
587 %!testif HAVE_CHOLMOD
588 %! C = gallery ("tridiag", 5);
589 %! sparse_chol2inv (C, eps*10);
590 
591 %!testif HAVE_CHOLMOD
592 %! D = gallery ("wathen", 1, 1);
593 %! sparse_chol2inv (D, eps*10^4);
594 
595 */
596 
597 DEFUN (cholupdate, args, nargout,
598  doc: /* -*- texinfo -*-
599 @deftypefn {} {[@var{R1}, @var{info}] =} cholupdate (@var{R}, @var{u}, @var{op})
600 Update or downdate a Cholesky@tie{}factorization.
601 
602 Given an upper triangular matrix @var{R} and a column vector @var{u},
603 attempt to determine another upper triangular matrix @var{R1} such that
604 
605 @itemize @bullet
606 @item
607 @var{R1}'*@var{R1} = @var{R}'*@var{R} + @var{u}*@var{u}'
608 if @var{op} is @qcode{"+"}
609 
610 @item
611 @var{R1}'*@var{R1} = @var{R}'*@var{R} - @var{u}*@var{u}'
612 if @var{op} is @qcode{"-"}
613 @end itemize
614 
615 If @var{op} is @qcode{"-"}, @var{info} is set to
616 
617 @itemize
618 @item 0 if the downdate was successful,
619 
620 @item 1 if @var{R}'*@var{R} - @var{u}*@var{u}' is not positive definite,
621 
622 @item 2 if @var{R} is singular.
623 @end itemize
624 
625 If @var{info} is not present, an error message is printed in cases 1 and 2.
626 @seealso{chol, cholinsert, choldelete, cholshift}
627 @end deftypefn */)
628 {
629  int nargin = args.length ();
630 
631  if (nargin < 2 || nargin > 3)
632  print_usage ();
633 
634  octave_value argr = args(0);
635  octave_value argu = args(1);
636 
637  if (! argr.isnumeric () || ! argu.isnumeric ()
638  || (nargin > 2 && ! args(2).is_string ()))
639  print_usage ();
640 
641  octave_value_list retval (nargout == 2 ? 2 : 1);
642 
643  octave_idx_type n = argr.rows ();
644 
645  std::string op = (nargin < 3) ? "+" : args(2).string_value ();
646 
647  bool down = (op == "-");
648 
649  if (! down && op != "+")
650  error (R"(cholupdate: OP must be "+" or "-")");
651 
652  if (argr.columns () != n || argu.rows () != n || argu.columns () != 1)
653  error ("cholupdate: dimension mismatch between R and U");
654 
655  int err = 0;
656  if (argr.is_single_type () || argu.is_single_type ())
657  {
658  if (argr.isreal () && argu.isreal ())
659  {
660  // real case
661  FloatMatrix R = argr.float_matrix_value ();
663 
664  math::chol<FloatMatrix> fact;
665  fact.set (R);
666 
667  if (down)
668  err = fact.downdate (u);
669  else
670  fact.update (u);
671 
672  retval = ovl (get_chol_r (fact));
673  }
674  else
675  {
676  // complex case
680 
681  math::chol<FloatComplexMatrix> fact;
682  fact.set (R);
683 
684  if (down)
685  err = fact.downdate (u);
686  else
687  fact.update (u);
688 
689  retval = ovl (get_chol_r (fact));
690  }
691  }
692  else
693  {
694  if (argr.isreal () && argu.isreal ())
695  {
696  // real case
697  Matrix R = argr.matrix_value ();
698  ColumnVector u = argu.column_vector_value ();
699 
700  math::chol<Matrix> fact;
701  fact.set (R);
702 
703  if (down)
704  err = fact.downdate (u);
705  else
706  fact.update (u);
707 
708  retval = ovl (get_chol_r (fact));
709  }
710  else
711  {
712  // complex case
715 
716  math::chol<ComplexMatrix> fact;
717  fact.set (R);
718 
719  if (down)
720  err = fact.downdate (u);
721  else
722  fact.update (u);
723 
724  retval = ovl (get_chol_r (fact));
725  }
726  }
727 
728  if (nargout > 1)
729  retval(1) = err;
730  else if (err == 1)
731  error ("cholupdate: downdate violates positiveness");
732  else if (err == 2)
733  error ("cholupdate: singular matrix");
734 
735  return retval;
736 }
737 
738 /*
739 %!shared A, u, Ac, uc
740 %! A = [ 0.436997 -0.131721 0.124120 -0.061673 ;
741 %! -0.131721 0.738529 0.019851 -0.140295 ;
742 %! 0.124120 0.019851 0.354879 -0.059472 ;
743 %! -0.061673 -0.140295 -0.059472 0.600939 ];
744 %!
745 %! u = [ 0.98950 ;
746 %! 0.39844 ;
747 %! 0.63484 ;
748 %! 0.13351 ];
749 %! Ac = [ 0.5585528 + 0.0000000i -0.1662088 - 0.0315341i 0.0107873 + 0.0236411i -0.0276775 - 0.0186073i ;
750 %! -0.1662088 + 0.0315341i 0.6760061 + 0.0000000i 0.0011452 - 0.0475528i 0.0145967 + 0.0247641i ;
751 %! 0.0107873 - 0.0236411i 0.0011452 + 0.0475528i 0.6263149 - 0.0000000i -0.1585837 - 0.0719763i ;
752 %! -0.0276775 + 0.0186073i 0.0145967 - 0.0247641i -0.1585837 + 0.0719763i 0.6034234 - 0.0000000i ];
753 %!
754 %! uc = [ 0.54267 + 0.91519i ;
755 %! 0.99647 + 0.43141i ;
756 %! 0.83760 + 0.68977i ;
757 %! 0.39160 + 0.90378i ];
758 
759 %!test
760 %! R = chol (A);
761 %! R1 = cholupdate (R, u);
762 %! assert (norm (triu (R1)-R1, Inf), 0);
763 %! assert (norm (R1'*R1 - R'*R - u*u', Inf), 0, 1e1*eps);
764 %!
765 %! R1 = cholupdate (R1, u, "-");
766 %! assert (norm (triu (R1)-R1, Inf), 0);
767 %! assert (norm (R1 - R, Inf), 0, 1e1*eps);
768 
769 %!test
770 %! R = chol (Ac);
771 %! R1 = cholupdate (R, uc);
772 %! assert (norm (triu (R1)-R1, Inf), 0);
773 %! assert (norm (R1'*R1 - R'*R - uc*uc', Inf), 0, 1e1*eps);
774 %!
775 %! R1 = cholupdate (R1, uc, "-");
776 %! assert (norm (triu (R1)-R1, Inf), 0);
777 %! assert (norm (R1 - R, Inf), 0, 11*eps);
778 
779 %!test
780 %! R = chol (single (A));
781 %! R1 = cholupdate (R, single (u));
782 %! assert (norm (triu (R1)-R1, Inf), single (0));
783 %! assert (norm (R1'*R1 - R'*R - single (u*u'), Inf), 0, 1e1* eps ("single"));
784 %!
785 %! R1 = cholupdate (R1, single (u), "-");
786 %! assert (norm (triu (R1)-R1, Inf), single (0));
787 %! assert (norm (R1 - R, Inf), 0, 2e1* eps ("single"));
788 
789 %!test
790 %! R = chol (single (Ac));
791 %! R1 = cholupdate (R, single (uc));
792 %! assert (norm (triu (R1)-R1, Inf), single (0));
793 %! assert (norm (R1'*R1 - R'*R - single (uc*uc'), Inf), 0, 1e1* eps ("single"));
794 %!
795 %! R1 = cholupdate (R1, single (uc), "-");
796 %! assert (norm (triu (R1)-R1, Inf), single (0));
797 %! assert (norm (R1 - R, Inf), 0, 2e1* eps ("single"));
798 */
799 
800 DEFUN (cholinsert, args, nargout,
801  doc: /* -*- texinfo -*-
802 @deftypefn {} {@var{R1} =} cholinsert (@var{R}, @var{j}, @var{u})
803 @deftypefnx {} {[@var{R1}, @var{info}] =} cholinsert (@var{R}, @var{j}, @var{u})
804 Update a Cholesky factorization given a row or column to insert in the
805 original factored matrix.
806 
807 Given a Cholesky@tie{}factorization of a real symmetric or complex Hermitian
808 positive definite matrix @w{@var{A} = @var{R}'*@var{R}}, @var{R}@tie{}upper
809 triangular, return the Cholesky@tie{}factorization of
810 @var{A1}, where @w{A1(p,p) = A}, @w{A1(:,j) = A1(j,:)' = u} and
811 @w{p = [1:j-1,j+1:n+1]}. @w{u(j)} should be positive.
812 
813 On return, @var{info} is set to
814 
815 @itemize
816 @item 0 if the insertion was successful,
817 
818 @item 1 if @var{A1} is not positive definite,
819 
820 @item 2 if @var{R} is singular.
821 @end itemize
822 
823 If @var{info} is not present, an error message is printed in cases 1 and 2.
824 @seealso{chol, cholupdate, choldelete, cholshift}
825 @end deftypefn */)
826 {
827  if (args.length () != 3)
828  print_usage ();
829 
830  octave_value argr = args(0);
831  octave_value argj = args(1);
832  octave_value argu = args(2);
833 
834  if (! argr.isnumeric () || ! argu.isnumeric ()
835  || ! argj.is_real_scalar ())
836  print_usage ();
837 
838  octave_idx_type n = argr.rows ();
839  octave_idx_type j = argj.scalar_value ();
840 
841  if (argr.columns () != n || argu.rows () != n+1 || argu.columns () != 1)
842  error ("cholinsert: dimension mismatch between R and U");
843 
844  if (j < 1 || j > n+1)
845  error ("cholinsert: index J out of range");
846 
847  octave_value_list retval (nargout == 2 ? 2 : 1);
848 
849  int err = 0;
850  if (argr.is_single_type () || argu.is_single_type ())
851  {
852  if (argr.isreal () && argu.isreal ())
853  {
854  // real case
855  FloatMatrix R = argr.float_matrix_value ();
857 
858  math::chol<FloatMatrix> fact;
859  fact.set (R);
860  err = fact.insert_sym (u, j-1);
861 
862  retval = ovl (get_chol_r (fact));
863  }
864  else
865  {
866  // complex case
870 
871  math::chol<FloatComplexMatrix> fact;
872  fact.set (R);
873  err = fact.insert_sym (u, j-1);
874 
875  retval = ovl (get_chol_r (fact));
876  }
877  }
878  else
879  {
880  if (argr.isreal () && argu.isreal ())
881  {
882  // real case
883  Matrix R = argr.matrix_value ();
884  ColumnVector u = argu.column_vector_value ();
885 
886  math::chol<Matrix> fact;
887  fact.set (R);
888  err = fact.insert_sym (u, j-1);
889 
890  retval = ovl (get_chol_r (fact));
891  }
892  else
893  {
894  // complex case
897 
898  math::chol<ComplexMatrix> fact;
899  fact.set (R);
900  err = fact.insert_sym (u, j-1);
901 
902  retval = ovl (get_chol_r (fact));
903  }
904  }
905 
906  if (nargout > 1)
907  retval(1) = err;
908  else if (err == 1)
909  error ("cholinsert: insertion violates positiveness");
910  else if (err == 2)
911  error ("cholinsert: singular matrix");
912  else if (err == 3)
913  error ("cholinsert: diagonal element must be real");
914 
915  return retval;
916 }
917 
918 /*
919 %!test
920 %! u2 = [ 0.35080 ;
921 %! 0.63930 ;
922 %! 3.31057 ;
923 %! -0.13825 ;
924 %! 0.45266 ];
925 %!
926 %! R = chol (A);
927 %!
928 %! j = 3; p = [1:j-1, j+1:5];
929 %! R1 = cholinsert (R, j, u2);
930 %! A1 = R1'*R1;
931 %!
932 %! assert (norm (triu (R1)-R1, Inf), 0);
933 %! assert (norm (A1(p,p) - A, Inf), 0, 1e1*eps);
934 
935 %!test
936 %! u2 = [ 0.35080 + 0.04298i;
937 %! 0.63930 + 0.23778i;
938 %! 3.31057 + 0.00000i;
939 %! -0.13825 + 0.19879i;
940 %! 0.45266 + 0.50020i];
941 %!
942 %! R = chol (Ac);
943 %!
944 %! j = 3; p = [1:j-1, j+1:5];
945 %! R1 = cholinsert (R, j, u2);
946 %! A1 = R1'*R1;
947 %!
948 %! assert (norm (triu (R1)-R1, Inf), 0);
949 %! assert (norm (A1(p,p) - Ac, Inf), 0, 1e1*eps);
950 
951 %!test
952 %! u2 = single ([ 0.35080 ;
953 %! 0.63930 ;
954 %! 3.31057 ;
955 %! -0.13825 ;
956 %! 0.45266 ]);
957 %!
958 %! R = chol (single (A));
959 %!
960 %! j = 3; p = [1:j-1, j+1:5];
961 %! R1 = cholinsert (R, j, u2);
962 %! A1 = R1'*R1;
963 %!
964 %! assert (norm (triu (R1)-R1, Inf), single (0));
965 %! assert (norm (A1(p,p) - A, Inf), 0, 1e1* eps ("single"));
966 
967 %!test
968 %! u2 = single ([ 0.35080 + 0.04298i;
969 %! 0.63930 + 0.23778i;
970 %! 3.31057 + 0.00000i;
971 %! -0.13825 + 0.19879i;
972 %! 0.45266 + 0.50020i]);
973 %!
974 %! R = chol (single (Ac));
975 %!
976 %! j = 3; p = [1:j-1, j+1:5];
977 %! R1 = cholinsert (R, j, u2);
978 %! A1 = R1'*R1;
979 %!
980 %! assert (norm (triu (R1)-R1, Inf), single (0));
981 %! assert (norm (A1(p,p) - single (Ac), Inf), 0, 2e1* eps ("single"));
982 
983 %!test
984 %! cu = chol (triu (A), "upper");
985 %! cl = chol (tril (A), "lower");
986 %! assert (cu, cl', eps);
987 
988 %!test
989 %! cca = chol (Ac);
990 %!
991 %! ccal = chol (Ac, "lower");
992 %! ccal2 = chol (tril (Ac), "lower");
993 %!
994 %! ccau = chol (Ac, "upper");
995 %! ccau2 = chol (triu (Ac), "upper");
996 %!
997 %! assert (cca'*cca, Ac, eps);
998 %! assert (ccau'*ccau, Ac, eps);
999 %! assert (ccau2'*ccau2, Ac, eps);
1000 %!
1001 %! assert (cca, ccal', eps);
1002 %! assert (cca, ccau, eps);
1003 %! assert (cca, ccal2', eps);
1004 %! assert (cca, ccau2, eps);
1005 
1006 %!test
1007 %! cca = chol (single (Ac));
1008 %!
1009 %! ccal = chol (single (Ac), "lower");
1010 %! ccal2 = chol (tril (single (Ac)), "lower");
1011 %!
1012 %! ccau = chol (single (Ac), "upper");
1013 %! ccau2 = chol (triu (single (Ac)), "upper");
1014 %!
1015 %! assert (cca'*cca, single (Ac), eps ("single"));
1016 %! assert (ccau'*ccau, single (Ac), eps ("single"));
1017 %! assert (ccau2'*ccau2, single (Ac), eps ("single"));
1018 %!
1019 %! assert (cca, ccal', eps ("single"));
1020 %! assert (cca, ccau, eps ("single"));
1021 %! assert (cca, ccal2', eps ("single"));
1022 %! assert (cca, ccau2, eps ("single"));
1023 
1024 %!test
1025 %! a = [12, 2, 3, 4;
1026 %! 2, 14, 5, 3;
1027 %! 3, 5, 16, 6;
1028 %! 4, 3, 6, 16];
1029 %!
1030 %! b = [0, 1, 2, 3;
1031 %! -1, 0, 1, 2;
1032 %! -2, -1, 0, 1;
1033 %! -3, -2, -1, 0];
1034 %!
1035 %! ca = a + i*b;
1036 %!
1037 %! cca = chol (ca);
1038 %!
1039 %! ccal = chol (ca, "lower");
1040 %! ccal2 = chol (tril (ca), "lower");
1041 %!
1042 %! ccau = chol (ca, "upper");
1043 %! ccau2 = chol (triu (ca), "upper");
1044 %!
1045 %! assert (cca'*cca, ca, 16*eps);
1046 %! assert (ccau'*ccau, ca, 16*eps);
1047 %! assert (ccau2'*ccau2, ca, 16*eps);
1048 %!
1049 %! assert (cca, ccal', 16*eps);
1050 %! assert (cca, ccau, 16*eps);
1051 %! assert (cca, ccal2', 16*eps);
1052 %! assert (cca, ccau2, 16*eps);
1053 */
1054 
1055 DEFUN (choldelete, args, ,
1056  doc: /* -*- texinfo -*-
1057 @deftypefn {} {@var{R1} =} choldelete (@var{R}, @var{j})
1058 Update a Cholesky factorization given a row or column to delete from the
1059 original factored matrix.
1060 
1061 Given a Cholesky@tie{}factorization of a real symmetric or complex Hermitian
1062 positive definite matrix @w{@var{A} = @var{R}'*@var{R}}, @var{R}@tie{}upper
1063 triangular, return the Cholesky@tie{}factorization of @w{A(p,p)}, where
1064 @w{p = [1:j-1,j+1:n+1]}.
1065 @seealso{chol, cholupdate, cholinsert, cholshift}
1066 @end deftypefn */)
1067 {
1068  if (args.length () != 2)
1069  print_usage ();
1070 
1071  octave_value argr = args(0);
1072  octave_value argj = args(1);
1073 
1074  if (! argr.isnumeric () || ! argj.is_real_scalar ())
1075  print_usage ();
1076 
1077  octave_idx_type n = argr.rows ();
1078  octave_idx_type j = argj.scalar_value ();
1079 
1080  if (argr.columns () != n)
1081  err_square_matrix_required ("choldelete", "R");
1082 
1083  if (j < 0 && j > n)
1084  error ("choldelete: index J out of range");
1085 
1086  octave_value_list retval;
1087 
1088  if (argr.is_single_type ())
1089  {
1090  if (argr.isreal ())
1091  {
1092  // real case
1093  FloatMatrix R = argr.float_matrix_value ();
1094 
1095  math::chol<FloatMatrix> fact;
1096  fact.set (R);
1097  fact.delete_sym (j-1);
1098 
1099  retval = ovl (get_chol_r (fact));
1100  }
1101  else
1102  {
1103  // complex case
1105 
1106  math::chol<FloatComplexMatrix> fact;
1107  fact.set (R);
1108  fact.delete_sym (j-1);
1109 
1110  retval = ovl (get_chol_r (fact));
1111  }
1112  }
1113  else
1114  {
1115  if (argr.isreal ())
1116  {
1117  // real case
1118  Matrix R = argr.matrix_value ();
1119 
1120  math::chol<Matrix> fact;
1121  fact.set (R);
1122  fact.delete_sym (j-1);
1123 
1124  retval = ovl (get_chol_r (fact));
1125  }
1126  else
1127  {
1128  // complex case
1129  ComplexMatrix R = argr.complex_matrix_value ();
1130 
1131  math::chol<ComplexMatrix> fact;
1132  fact.set (R);
1133  fact.delete_sym (j-1);
1134 
1135  retval = ovl (get_chol_r (fact));
1136  }
1137  }
1138 
1139  return retval;
1140 }
1141 
1142 /*
1143 %!test
1144 %! R = chol (A);
1145 %!
1146 %! j = 3; p = [1:j-1,j+1:4];
1147 %! R1 = choldelete (R, j);
1148 %!
1149 %! assert (norm (triu (R1)-R1, Inf), 0);
1150 %! assert (norm (R1'*R1 - A(p,p), Inf), 0, 1e1*eps);
1151 
1152 %!test
1153 %! R = chol (Ac);
1154 %!
1155 %! j = 3; p = [1:j-1,j+1:4];
1156 %! R1 = choldelete (R, j);
1157 %!
1158 %! assert (norm (triu (R1)-R1, Inf), 0);
1159 %! assert (norm (R1'*R1 - Ac(p,p), Inf), 0, 1e1*eps);
1160 
1161 %!test
1162 %! R = chol (single (A));
1163 %!
1164 %! j = 3; p = [1:j-1,j+1:4];
1165 %! R1 = choldelete (R, j);
1166 %!
1167 %! assert (norm (triu (R1)-R1, Inf), single (0));
1168 %! assert (norm (R1'*R1 - single (A(p,p)), Inf), 0, 1e1* eps ("single"));
1169 
1170 %!test
1171 %! R = chol (single (Ac));
1172 %!
1173 %! j = 3; p = [1:j-1,j+1:4];
1174 %! R1 = choldelete (R,j);
1175 %!
1176 %! assert (norm (triu (R1)-R1, Inf), single (0));
1177 %! assert (norm (R1'*R1 - single (Ac(p,p)), Inf), 0, 1e1* eps ("single"));
1178 */
1179 
1180 DEFUN (cholshift, args, ,
1181  doc: /* -*- texinfo -*-
1182 @deftypefn {} {@var{R1} =} cholshift (@var{R}, @var{i}, @var{j})
1183 Update a Cholesky factorization given a range of columns to shift in the
1184 original factored matrix.
1185 
1186 Given a Cholesky@tie{}factorization of a real symmetric or complex Hermitian
1187 positive definite matrix @w{@var{A} = @var{R}'*@var{R}}, @var{R}@tie{}upper
1188 triangular, return the Cholesky@tie{}factorization of
1189 @w{@var{A}(p,p)}, where @w{p} is the permutation @*
1190 @code{p = [1:i-1, shift(i:j, 1), j+1:n]} if @w{@var{i} < @var{j}} @*
1191  or @*
1192 @code{p = [1:j-1, shift(j:i,-1), i+1:n]} if @w{@var{j} < @var{i}}. @*
1193 
1194 @seealso{chol, cholupdate, cholinsert, choldelete}
1195 @end deftypefn */)
1196 {
1197  if (args.length () != 3)
1198  print_usage ();
1199 
1200  octave_value argr = args(0);
1201  octave_value argi = args(1);
1202  octave_value argj = args(2);
1203 
1204  if (! argr.isnumeric () || ! argi.is_real_scalar ()
1205  || ! argj.is_real_scalar ())
1206  print_usage ();
1207 
1208  octave_idx_type n = argr.rows ();
1209  octave_idx_type i = argi.scalar_value ();
1210  octave_idx_type j = argj.scalar_value ();
1211 
1212  if (argr.columns () != n)
1213  err_square_matrix_required ("cholshift", "R");
1214 
1215  if (j < 0 || j > n+1 || i < 0 || i > n+1)
1216  error ("cholshift: index I or J is out of range");
1217 
1218  octave_value_list retval;
1219 
1220  if (argr.is_single_type () && argi.is_single_type ()
1221  && argj.is_single_type ())
1222  {
1223  if (argr.isreal ())
1224  {
1225  // real case
1226  FloatMatrix R = argr.float_matrix_value ();
1227 
1228  math::chol<FloatMatrix> fact;
1229  fact.set (R);
1230  fact.shift_sym (i-1, j-1);
1231 
1232  retval = ovl (get_chol_r (fact));
1233  }
1234  else
1235  {
1236  // complex case
1238 
1239  math::chol<FloatComplexMatrix> fact;
1240  fact.set (R);
1241  fact.shift_sym (i-1, j-1);
1242 
1243  retval = ovl (get_chol_r (fact));
1244  }
1245  }
1246  else
1247  {
1248  if (argr.isreal ())
1249  {
1250  // real case
1251  Matrix R = argr.matrix_value ();
1252 
1253  math::chol<Matrix> fact;
1254  fact.set (R);
1255  fact.shift_sym (i-1, j-1);
1256 
1257  retval = ovl (get_chol_r (fact));
1258  }
1259  else
1260  {
1261  // complex case
1262  ComplexMatrix R = argr.complex_matrix_value ();
1263 
1264  math::chol<ComplexMatrix> fact;
1265  fact.set (R);
1266  fact.shift_sym (i-1, j-1);
1267 
1268  retval = ovl (get_chol_r (fact));
1269  }
1270  }
1271 
1272  return retval;
1273 }
1274 
1275 /*
1276 %!test
1277 %! R = chol (A);
1278 %!
1279 %! i = 1; j = 3; p = [1:i-1, shift(i:j,-1), j+1:4];
1280 %! R1 = cholshift (R, i, j);
1281 %!
1282 %! assert (norm (triu (R1)-R1, Inf), 0);
1283 %! assert (norm (R1'*R1 - A(p,p), Inf), 0, 1e1*eps);
1284 %!
1285 %! j = 1; i = 3; p = [1:j-1, shift(j:i,+1), i+1:4];
1286 %! R1 = cholshift (R, i, j);
1287 %!
1288 %! assert (norm (triu (R1) - R1, Inf), 0);
1289 %! assert (norm (R1'*R1 - A(p,p), Inf), 0, 1e1*eps);
1290 
1291 %!test
1292 %! R = chol (Ac);
1293 %!
1294 %! i = 1; j = 3; p = [1:i-1, shift(i:j,-1), j+1:4];
1295 %! R1 = cholshift (R, i, j);
1296 %!
1297 %! assert (norm (triu (R1)-R1, Inf), 0);
1298 %! assert (norm (R1'*R1 - Ac(p,p), Inf), 0, 1e1*eps);
1299 %!
1300 %! j = 1; i = 3; p = [1:j-1, shift(j:i,+1), i+1:4];
1301 %! R1 = cholshift (R, i, j);
1302 %!
1303 %! assert (norm (triu (R1)-R1, Inf), 0);
1304 %! assert (norm (R1'*R1 - Ac(p,p), Inf), 0, 1e1*eps);
1305 
1306 %!test
1307 %! R = chol (single (A));
1308 %!
1309 %! i = 1; j = 3; p = [1:i-1, shift(i:j,-1), j+1:4];
1310 %! R1 = cholshift (R, i, j);
1311 %!
1312 %! assert (norm (triu (R1)-R1, Inf), 0);
1313 %! assert (norm (R1'*R1 - single (A(p,p)), Inf), 0, 1e1* eps ("single"));
1314 %!
1315 %! j = 1; i = 3; p = [1:j-1, shift(j:i,+1), i+1:4];
1316 %! R1 = cholshift (R, i, j);
1317 %!
1318 %! assert (norm (triu (R1)-R1, Inf), 0);
1319 %! assert (norm (R1'*R1 - single (A(p,p)), Inf), 0, 1e1* eps ("single"));
1320 
1321 %!test
1322 %! R = chol (single (Ac));
1323 %!
1324 %! i = 1; j = 3; p = [1:i-1, shift(i:j,-1), j+1:4];
1325 %! R1 = cholshift (R, i, j);
1326 %!
1327 %! assert (norm (triu (R1)-R1, Inf), 0);
1328 %! assert (norm (R1'*R1 - single (Ac(p,p)), Inf), 0, 1e1* eps ("single"));
1329 %!
1330 %! j = 1; i = 3; p = [1:j-1, shift(j:i,+1), i+1:4];
1331 %! R1 = cholshift (R, i, j);
1332 %!
1333 %! assert (norm (triu (R1)-R1, Inf), 0);
1334 %! assert (norm (R1'*R1 - single (Ac(p,p)), Inf), 0, 1e1* eps ("single"));
1335 */
1336 
1337 OCTAVE_END_NAMESPACE(octave)
Definition: dMatrix.h:42
Definition: chol.h:38
T inverse() const
Definition: chol.cc:250
SparseMatrix sparse_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:900
bool is_real_scalar() const
Definition: ov.h:610
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_single_type() const
Definition: ov.h:698
bool isempty() const
Definition: ov.h:601
FloatColumnVector float_column_vector_value(bool frc_str_conv=false, bool frc_vec_conv=false) const
bool issparse() const
Definition: ov.h:753
double scalar_value(bool frc_str_conv=false) const
Definition: ov.h:847
FloatComplexColumnVector float_complex_column_vector_value(bool frc_str_conv=false, bool frc_vec_conv=false) const
FloatMatrix float_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:856
bool iscomplex() const
Definition: ov.h:741
ComplexColumnVector complex_column_vector_value(bool frc_str_conv=false, bool frc_vec_conv=false) const
ColumnVector column_vector_value(bool frc_str_conv=false, bool frc_vec_conv=false) const
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
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() error(const char *fmt,...)
Definition: error.cc:988
void err_square_matrix_required(const char *fcn, const char *name)
Definition: errwarn.cc:122
void err_wrong_type_arg(const char *name, const char *s)
Definition: errwarn.cc:166
T chol2inv(const T &r)
Definition: chol.cc:242
T octave_idx_type m
Definition: mx-inlines.cc:781
octave_idx_type n
Definition: mx-inlines.cc:761
T * r
Definition: mx-inlines.cc:781
bool strcmpi(const T &str_a, const T &str_b)
True if strings are the same, ignoring case.
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