47 template <
class CHOLT>
55 template <
class CHOLT>
65 @deftypefn {Loadable Function} {@var{R} =} chol (@var{A})\n\
66 @deftypefnx {Loadable Function} {[@var{R}, @var{p}] =} chol (@var{A})\n\
67 @deftypefnx {Loadable Function} {[@var{R}, @var{p}, @var{Q}] =} chol (@var{S})\n\
68 @deftypefnx {Loadable Function} {[@var{R}, @var{p}, @var{Q}] =} chol (@var{S}, \"vector\")\n\
69 @deftypefnx {Loadable Function} {[@var{L}, @dots{}] =} chol (@dots{}, \"lower\")\n\
70 @deftypefnx {Loadable Function} {[@var{L}, @dots{}] =} chol (@dots{}, \"upper\")\n\
71 @cindex Cholesky factorization\n\
72 Compute the Cholesky@tie{}factor, @var{R}, of the symmetric positive definite\n\
73 matrix @var{A}, where\n\
80 @var{R}' * @var{R} = @var{A}.\n\
85 Called with one output argument @code{chol} fails if @var{A} or @var{S} is\n\
86 not positive definite. With two or more output arguments @var{p} flags\n\
87 whether the matrix was positive definite and @code{chol} does not fail. A\n\
88 zero value indicated that the matrix was positive definite and the @var{R}\n\
89 gives the factorization, and @var{p} will have a positive value otherwise.\n\
91 If called with 3 outputs then a sparsity preserving row/column permutation\n\
92 is applied to @var{A} prior to the factorization. That is @var{R}\n\
93 is the factorization of @code{@var{A}(@var{Q},@var{Q})} such that\n\
95 $ R^T R = Q^T A Q$.\n\
100 @var{R}' * @var{R} = @var{Q}' * @var{A} * @var{Q}.\n\
105 The sparsity preserving permutation is generally returned as a matrix.\n\
106 However, given the flag @qcode{\"vector\"}, @var{Q} will be returned as a\n\
109 $ R^T R = A (Q, Q)$.\n\
114 @var{R}' * @var{R} = @var{A}(@var{Q}, @var{Q}).\n\
119 Called with either a sparse or full matrix and using the @qcode{\"lower\"}\n\
120 flag, @code{chol} returns the lower triangular factorization such that\n\
127 @var{L} * @var{L}' = @var{A}.\n\
132 For full matrices, if the @qcode{\"lower\"} flag is set only the lower\n\
133 triangular part of the matrix is used for the factorization, otherwise the\n\
134 upper triangular part is used.\n\
136 In general the lower triangular factorization is significantly faster for\n\
138 @seealso{hess, lu, qr, qz, schur, svd, cholinv, chol2inv, cholupdate, cholinsert, choldelete, cholshift}\n\
142 int nargin = args.
length ();
146 if (nargin < 1 || nargin > 3 || nargout > 3
147 || (! args(0).is_sparse_type () && nargout > 2))
156 std::string tmp = args(n++).string_value ();
160 if (tmp.compare (
"vector") == 0)
162 else if (tmp.compare (
"lower") == 0)
168 else if (tmp.compare (
"upper") == 0)
171 error (
"chol: unexpected second or third input");
174 error (
"chol: expecting trailing string arguments");
184 int arg_is_empty =
empty_arg (
"chol", nr, nc);
186 if (arg_is_empty < 0)
188 if (arg_is_empty > 0)
194 bool natural = (nargout != 3);
195 bool force = nargout > 1;
208 retval(2) = fact.
perm ();
210 retval(2) = fact.
Q ();
213 if (nargout > 1 || info == 0)
215 retval(1) = fact.
P ();
217 retval(0) = fact.
L ();
219 retval(0) = fact.
R ();
222 error (
"chol: input matrix must be positive definite");
236 retval(2) = fact.
perm ();
238 retval(2) = fact.
Q ();
241 if (nargout > 1 || info == 0)
243 retval(1) = fact.
P ();
245 retval(0) = fact.
L ();
247 retval(0) = fact.
R ();
250 error (
"chol: input matrix must be positive definite");
272 if (nargout == 2 || info == 0)
281 error (
"chol: input matrix must be positive definite");
298 if (nargout == 2 || info == 0)
307 error (
"chol: input matrix must be positive definite");
327 fact =
CHOL (m, info);
329 if (nargout == 2 || info == 0)
338 error (
"chol: input matrix must be positive definite");
355 if (nargout == 2 || info == 0)
364 error (
"chol: input matrix must be positive definite");
387 @deftypefn {Loadable Function} {} cholinv (@var{A})\n\
388 Use the Cholesky@tie{}factorization to compute the inverse of the\n\
389 symmetric positive definite matrix @var{A}.\n\
390 @seealso{chol, chol2inv, inv}\n\
395 int nargin = args.
length ();
404 if (nr == 0 || nc == 0)
423 error (
"cholinv: A must be positive definite");
437 error (
"cholinv: A must be positive definite");
456 error (
"cholinv: A must be positive definite");
470 error (
"cholinv: A must be positive definite");
489 error (
"cholinv: A must be positive definite");
503 error (
"cholinv: A must be positive definite");
534 @deftypefn {Loadable Function} {} chol2inv (@var{U})\n\
535 Invert a symmetric, positive definite square matrix from its Cholesky\n\
536 decomposition, @var{U}. Note that @var{U} should be an upper-triangular\n\
537 matrix with positive diagonal elements. @code{chol2inv (@var{U})}\n\
538 provides @code{inv (@var{U}'*@var{U})} but it is much faster than\n\
540 @seealso{chol, cholinv, inv}\n\
545 int nargin = args.
length ();
554 if (nr == 0 || nc == 0)
626 @deftypefn {Loadable Function} {[@var{R1}, @var{info}] =} cholupdate (@var{R}, @var{u}, @var{op})\n\
627 Update or downdate a Cholesky@tie{}factorization. Given an upper triangular\n\
628 matrix @var{R} and a column vector @var{u}, attempt to determine another\n\
629 upper triangular matrix @var{R1} such that\n\
633 @var{R1}'*@var{R1} = @var{R}'*@var{R} + @var{u}*@var{u}'\n\
634 if @var{op} is @qcode{\"+\"}\n\
637 @var{R1}'*@var{R1} = @var{R}'*@var{R} - @var{u}*@var{u}'\n\
638 if @var{op} is @qcode{\"-\"}\n\
641 If @var{op} is @qcode{\"-\"}, @var{info} is set to\n\
644 @item 0 if the downdate was successful,\n\
646 @item 1 if @var{R}'*@var{R} - @var{u}*@var{u}' is not positive definite,\n\
648 @item 2 if @var{R} is singular.\n\
651 If @var{info} is not present, an error message is printed in cases 1 and 2.\n\
652 @seealso{chol, cholinsert, choldelete, cholshift}\n\
659 if (nargin > 3 || nargin < 2)
669 && (nargin < 3 || args(2).is_string ()))
673 std::string op = (nargin < 3) ?
"+" : args(2).string_value ();
675 bool down = op ==
"-";
677 if (down || op ==
"+")
756 error (
"cholupdate: downdate violates positiveness");
758 error (
"cholupdate: singular matrix");
761 error (
"cholupdate: dimension mismatch between R and U");
763 error (
"cholupdate: OP must be \"+\" or \"-\"");
835 @deftypefn {Loadable Function} {@var{R1} =} cholinsert (@var{R}, @var{j}, @var{u})\n\
836 @deftypefnx {Loadable Function} {[@var{R1}, @var{info}] =} cholinsert (@var{R}, @var{j}, @var{u})\n\
837 Given a Cholesky@tie{}factorization of a real symmetric or complex Hermitian\n\
838 positive definite matrix @w{@var{A} = @var{R}'*@var{R}}, @var{R}@tie{}upper\n\
839 triangular, return the Cholesky@tie{}factorization of\n\
840 @var{A1}, where @w{A1(p,p) = A}, @w{A1(:,j) = A1(j,:)' = u} and\n\
841 @w{p = [1:j-1,j+1:n+1]}. @w{u(j)} should be positive.\n\
842 On return, @var{info} is set to\n\
845 @item 0 if the insertion was successful,\n\
847 @item 1 if @var{A1} is not positive definite,\n\
849 @item 2 if @var{R} is singular.\n\
852 If @var{info} is not present, an error message is printed in cases 1 and 2.\n\
853 @seealso{chol, cholupdate, choldelete, cholshift}\n\
878 if (j > 0 && j <= n+1)
941 error (
"cholinsert: insertion violates positiveness");
943 error (
"cholinsert: singular matrix");
945 error (
"cholinsert: diagonal element must be real");
948 error (
"cholinsert: index J out of range");
951 error (
"cholinsert: dimension mismatch between R and U");
1098 @deftypefn {Loadable Function} {@var{R1} =} choldelete (@var{R}, @var{j})\n\
1099 Given a Cholesky@tie{}factorization of a real symmetric or complex Hermitian\n\
1100 positive definite matrix @w{@var{A} = @var{R}'*@var{R}}, @var{R}@tie{}upper\n\
1101 triangular, return the Cholesky@tie{}factorization of @w{A(p,p)}, where\n\
1102 @w{p = [1:j-1,j+1:n+1]}.\n\
1103 @seealso{chol, cholupdate, cholinsert, cholshift}\n\
1126 if (j > 0 && j <= n)
1180 error (
"choldelete: index J out of range");
1183 error (
"choldelete: matrix R must be square");
1231 @deftypefn {Loadable Function} {@var{R1} =} cholshift (@var{R}, @var{i}, @var{j})\n\
1232 Given a Cholesky@tie{}factorization of a real symmetric or complex Hermitian\n\
1233 positive definite matrix @w{@var{A} = @var{R}'*@var{R}}, @var{R}@tie{}upper\n\
1234 triangular, return the Cholesky@tie{}factorization of\n\
1235 @w{@var{A}(p,p)}, where @w{p} is the permutation @*\n\
1236 @code{p = [1:i-1, shift(i:j, 1), j+1:n]} if @w{@var{i} < @var{j}} @*\n\
1238 @code{p = [1:j-1, shift(j:i,-1), i+1:n]} if @w{@var{j} < @var{i}}. @*\n\
1240 @seealso{chol, cholupdate, cholinsert, choldelete}\n\
1266 if (j > 0 && j <= n+1 && i > 0 && i <= n+1)
1322 error (
"cholshift: index I or J is out of range");
1325 error (
"cholshift: R must be a square matrix");