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