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