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