GNU Octave  3.8.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
chol.cc
Go to the documentation of this file.
1 /*
2 
3 Copyright (C) 1996-2013 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 the
11 Free Software Foundation; either version 3 of the License, or (at your
12 option) any later version.
13 
14 Octave is distributed in the hope that it will be useful, but WITHOUT
15 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
16 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
17 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 <http://www.gnu.org/licenses/>.
22 
23 */
24 
25 
26 #ifdef HAVE_CONFIG_H
27 #include <config.h>
28 #endif
29 
30 #include "CmplxCHOL.h"
31 #include "dbleCHOL.h"
32 #include "fCmplxCHOL.h"
33 #include "floatCHOL.h"
34 #include "SparseCmplxCHOL.h"
35 #include "SparsedbleCHOL.h"
36 #include "oct-spparms.h"
37 #include "sparse-util.h"
38 
39 #include "ov-re-sparse.h"
40 #include "ov-cx-sparse.h"
41 #include "defun-dld.h"
42 #include "error.h"
43 #include "gripes.h"
44 #include "oct-obj.h"
45 #include "utils.h"
46 
47 template <class CHOLT>
48 static octave_value
49 get_chol_r (const CHOLT& fact)
50 {
51  return octave_value (fact.chol_matrix (),
53 }
54 
55 template <class CHOLT>
56 static octave_value
57 get_chol_l (const CHOLT& fact)
58 {
59  return octave_value (fact.chol_matrix ().transpose (),
61 }
62 
63 DEFUN_DLD (chol, args, nargout,
64  "-*- texinfo -*-\n\
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\
74 @tex\n\
75 $ R^T R = A $.\n\
76 @end tex\n\
77 @ifnottex\n\
78 \n\
79 @example\n\
80 @var{R}' * @var{R} = @var{A}.\n\
81 @end example\n\
82 \n\
83 @end ifnottex\n\
84 \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\
90 \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\
94 @tex\n\
95 $ R^T R = Q^T A Q$.\n\
96 @end tex\n\
97 @ifnottex\n\
98 \n\
99 @example\n\
100 @var{R}' * @var{R} = @var{Q}' * @var{A} * @var{Q}.\n\
101 @end example\n\
102 \n\
103 @end ifnottex\n\
104 \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\
107 vector such that\n\
108 @tex\n\
109 $ R^T R = A (Q, Q)$.\n\
110 @end tex\n\
111 @ifnottex\n\
112 \n\
113 @example\n\
114 @var{R}' * @var{R} = @var{A}(@var{Q}, @var{Q}).\n\
115 @end example\n\
116 \n\
117 @end ifnottex\n\
118 \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\
121 @tex\n\
122 $ L L^T = A $.\n\
123 @end tex\n\
124 @ifnottex\n\
125 \n\
126 @example\n\
127 @var{L} * @var{L}' = @var{A}.\n\
128 @end example\n\
129 \n\
130 @end ifnottex\n\
131 \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\
135 \n\
136 In general the lower triangular factorization is significantly faster for\n\
137 sparse matrices.\n\
138 @seealso{hess, lu, qr, qz, schur, svd, cholinv, chol2inv, cholupdate, cholinsert, choldelete, cholshift}\n\
139 @end deftypefn")
140 {
141  octave_value_list retval;
142  int nargin = args.length ();
143  bool LLt = false;
144  bool vecout = false;
145 
146  if (nargin < 1 || nargin > 3 || nargout > 3
147  || (! args(0).is_sparse_type () && nargout > 2))
148  {
149  print_usage ();
150  return retval;
151  }
152 
153  int n = 1;
154  while (n < nargin && ! error_state)
155  {
156  std::string tmp = args(n++).string_value ();
157 
158  if (! error_state )
159  {
160  if (tmp.compare ("vector") == 0)
161  vecout = true;
162  else if (tmp.compare ("lower") == 0)
163  // FIXME: currently the option "lower" is handled by transposing
164  // the matrix, factorizing it with the lapack function
165  // DPOTRF ('U', ...) and finally transposing the factor. It would
166  // be more efficient to use DPOTRF ('L', ...) in this case.
167  LLt = true;
168  else if (tmp.compare ("upper") == 0)
169  LLt = false;
170  else
171  error ("chol: unexpected second or third input");
172  }
173  else
174  error ("chol: expecting trailing string arguments");
175  }
176 
177  if (! error_state)
178  {
179  octave_value arg = args(0);
180 
181  octave_idx_type nr = arg.rows ();
182  octave_idx_type nc = arg.columns ();
183 
184  int arg_is_empty = empty_arg ("chol", nr, nc);
185 
186  if (arg_is_empty < 0)
187  return retval;
188  if (arg_is_empty > 0)
189  return octave_value (Matrix ());
190 
191  if (arg.is_sparse_type ())
192  {
193  octave_idx_type info;
194  bool natural = (nargout != 3);
195  bool force = nargout > 1;
196 
197  if (arg.is_real_type ())
198  {
200 
201  if (! error_state)
202  {
203  SparseCHOL fact (m, info, natural, force);
204 
205  if (nargout == 3)
206  {
207  if (vecout)
208  retval(2) = fact.perm ();
209  else
210  retval(2) = fact.Q ();
211  }
212 
213  if (nargout > 1 || info == 0)
214  {
215  retval(1) = fact.P ();
216  if (LLt)
217  retval(0) = fact.L ();
218  else
219  retval(0) = fact.R ();
220  }
221  else
222  error ("chol: input matrix must be positive definite");
223  }
224  }
225  else if (arg.is_complex_type ())
226  {
228 
229  if (! error_state)
230  {
231  SparseComplexCHOL fact (m, info, natural, force);
232 
233  if (nargout == 3)
234  {
235  if (vecout)
236  retval(2) = fact.perm ();
237  else
238  retval(2) = fact.Q ();
239  }
240 
241  if (nargout > 1 || info == 0)
242  {
243  retval(1) = fact.P ();
244  if (LLt)
245  retval(0) = fact.L ();
246  else
247  retval(0) = fact.R ();
248  }
249  else
250  error ("chol: input matrix must be positive definite");
251  }
252  }
253  else
254  gripe_wrong_type_arg ("chol", arg);
255  }
256  else if (arg.is_single_type ())
257  {
258  if (arg.is_real_type ())
259  {
260  FloatMatrix m = arg.float_matrix_value ();
261 
262  if (! error_state)
263  {
264  octave_idx_type info;
265 
266  FloatCHOL fact;
267  if (LLt)
268  fact = FloatCHOL (m.transpose (), info);
269  else
270  fact = FloatCHOL (m, info);
271 
272  if (nargout == 2 || info == 0)
273  {
274  retval(1) = info;
275  if (LLt)
276  retval(0) = get_chol_l (fact);
277  else
278  retval(0) = get_chol_r (fact);
279  }
280  else
281  error ("chol: input matrix must be positive definite");
282  }
283  }
284  else if (arg.is_complex_type ())
285  {
287 
288  if (! error_state)
289  {
290  octave_idx_type info;
291 
292  FloatComplexCHOL fact;
293  if (LLt)
294  fact = FloatComplexCHOL (m.transpose (), info);
295  else
296  fact = FloatComplexCHOL (m, info);
297 
298  if (nargout == 2 || info == 0)
299  {
300  retval(1) = info;
301  if (LLt)
302  retval(0) = get_chol_l (fact);
303  else
304  retval(0) = get_chol_r (fact);
305  }
306  else
307  error ("chol: input matrix must be positive definite");
308  }
309  }
310  else
311  gripe_wrong_type_arg ("chol", arg);
312  }
313  else
314  {
315  if (arg.is_real_type ())
316  {
317  Matrix m = arg.matrix_value ();
318 
319  if (! error_state)
320  {
321  octave_idx_type info;
322 
323  CHOL fact;
324  if (LLt)
325  fact = CHOL (m.transpose (), info);
326  else
327  fact = CHOL (m, info);
328 
329  if (nargout == 2 || info == 0)
330  {
331  retval(1) = info;
332  if (LLt)
333  retval(0) = get_chol_l (fact);
334  else
335  retval(0) = get_chol_r (fact);
336  }
337  else
338  error ("chol: input matrix must be positive definite");
339  }
340  }
341  else if (arg.is_complex_type ())
342  {
344 
345  if (! error_state)
346  {
347  octave_idx_type info;
348 
349  ComplexCHOL fact;
350  if (LLt)
351  fact = ComplexCHOL (m.transpose (), info);
352  else
353  fact = ComplexCHOL (m, info);
354 
355  if (nargout == 2 || info == 0)
356  {
357  retval(1) = info;
358  if (LLt)
359  retval(0) = get_chol_l (fact);
360  else
361  retval(0) = get_chol_r (fact);
362  }
363  else
364  error ("chol: input matrix must be positive definite");
365  }
366  }
367  else
368  gripe_wrong_type_arg ("chol", arg);
369  }
370  }
371 
372  return retval;
373 }
374 
375 /*
376 %!assert (chol ([2, 1; 1, 1]), [sqrt(2), 1/sqrt(2); 0, 1/sqrt(2)], sqrt (eps))
377 %!assert (chol (single ([2, 1; 1, 1])), single ([sqrt(2), 1/sqrt(2); 0, 1/sqrt(2)]), sqrt (eps ("single")))
378 
379 %!error chol ()
380 %!error <matrix must be positive definite> chol ([1, 2; 3, 4])
381 %!error <requires square matrix> chol ([1, 2; 3, 4; 5, 6])
382 %!error <unexpected second or third input> chol (1, 2)
383 */
384 
385 DEFUN_DLD (cholinv, args, ,
386  "-*- texinfo -*-\n\
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\
391 @end deftypefn")
392 {
393  octave_value retval;
394 
395  int nargin = args.length ();
396 
397  if (nargin == 1)
398  {
399  octave_value arg = args(0);
400 
401  octave_idx_type nr = arg.rows ();
402  octave_idx_type nc = arg.columns ();
403 
404  if (nr == 0 || nc == 0)
405  retval = Matrix ();
406  else
407  {
408  if (arg.is_sparse_type ())
409  {
410  octave_idx_type info;
411 
412  if (arg.is_real_type ())
413  {
415 
416  if (! error_state)
417  {
418  SparseCHOL chol (m, info);
419 
420  if (info == 0)
421  retval = chol.inverse ();
422  else
423  error ("cholinv: A must be positive definite");
424  }
425  }
426  else if (arg.is_complex_type ())
427  {
429 
430  if (! error_state)
431  {
432  SparseComplexCHOL chol (m, info);
433 
434  if (info == 0)
435  retval = chol.inverse ();
436  else
437  error ("cholinv: A must be positive definite");
438  }
439  }
440  else
441  gripe_wrong_type_arg ("cholinv", arg);
442  }
443  else if (arg.is_single_type ())
444  {
445  if (arg.is_real_type ())
446  {
447  FloatMatrix m = arg.float_matrix_value ();
448 
449  if (! error_state)
450  {
451  octave_idx_type info;
452  FloatCHOL chol (m, info);
453  if (info == 0)
454  retval = chol.inverse ();
455  else
456  error ("cholinv: A must be positive definite");
457  }
458  }
459  else if (arg.is_complex_type ())
460  {
462 
463  if (! error_state)
464  {
465  octave_idx_type info;
466  FloatComplexCHOL chol (m, info);
467  if (info == 0)
468  retval = chol.inverse ();
469  else
470  error ("cholinv: A must be positive definite");
471  }
472  }
473  else
474  gripe_wrong_type_arg ("chol", arg);
475  }
476  else
477  {
478  if (arg.is_real_type ())
479  {
480  Matrix m = arg.matrix_value ();
481 
482  if (! error_state)
483  {
484  octave_idx_type info;
485  CHOL chol (m, info);
486  if (info == 0)
487  retval = chol.inverse ();
488  else
489  error ("cholinv: A must be positive definite");
490  }
491  }
492  else if (arg.is_complex_type ())
493  {
495 
496  if (! error_state)
497  {
498  octave_idx_type info;
499  ComplexCHOL chol (m, info);
500  if (info == 0)
501  retval = chol.inverse ();
502  else
503  error ("cholinv: A must be positive definite");
504  }
505  }
506  else
507  gripe_wrong_type_arg ("chol", arg);
508  }
509  }
510  }
511  else
512  print_usage ();
513 
514  return retval;
515 }
516 
517 /*
518 %!shared A, Ainv
519 %! A = [2,0.2;0.2,1];
520 %! Ainv = inv (A);
521 %!test
522 %! Ainv1 = cholinv (A);
523 %! assert (norm (Ainv-Ainv1), 0, 1e-10);
524 %!testif HAVE_CHOLMOD
525 %! Ainv2 = inv (sparse (A));
526 %! assert (norm (Ainv-Ainv2), 0, 1e-10);
527 %!testif HAVE_CHOLMOD
528 %! Ainv3 = cholinv (sparse (A));
529 %! assert (norm (Ainv-Ainv3), 0, 1e-10);
530 */
531 
533  "-*- texinfo -*-\n\
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\
539 using @code{inv}.\n\
540 @seealso{chol, cholinv, inv}\n\
541 @end deftypefn")
542 {
543  octave_value retval;
544 
545  int nargin = args.length ();
546 
547  if (nargin == 1)
548  {
549  octave_value arg = args(0);
550 
551  octave_idx_type nr = arg.rows ();
552  octave_idx_type nc = arg.columns ();
553 
554  if (nr == 0 || nc == 0)
555  retval = Matrix ();
556  else
557  {
558  if (arg.is_sparse_type ())
559  {
560  if (arg.is_real_type ())
561  {
563 
564  if (! error_state)
565  retval = chol2inv (r);
566  }
567  else if (arg.is_complex_type ())
568  {
570 
571  if (! error_state)
572  retval = chol2inv (r);
573  }
574  else
575  gripe_wrong_type_arg ("chol2inv", arg);
576  }
577  else if (arg.is_single_type ())
578  {
579  if (arg.is_real_type ())
580  {
581  FloatMatrix r = arg.float_matrix_value ();
582 
583  if (! error_state)
584  retval = chol2inv (r);
585  }
586  else if (arg.is_complex_type ())
587  {
589 
590  if (! error_state)
591  retval = chol2inv (r);
592  }
593  else
594  gripe_wrong_type_arg ("chol2inv", arg);
595 
596  }
597  else
598  {
599  if (arg.is_real_type ())
600  {
601  Matrix r = arg.matrix_value ();
602 
603  if (! error_state)
604  retval = chol2inv (r);
605  }
606  else if (arg.is_complex_type ())
607  {
609 
610  if (! error_state)
611  retval = chol2inv (r);
612  }
613  else
614  gripe_wrong_type_arg ("chol2inv", arg);
615  }
616  }
617  }
618  else
619  print_usage ();
620 
621  return retval;
622 }
623 
624 DEFUN_DLD (cholupdate, args, nargout,
625  "-*- texinfo -*-\n\
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\
630 \n\
631 @itemize @bullet\n\
632 @item\n\
633 @var{R1}'*@var{R1} = @var{R}'*@var{R} + @var{u}*@var{u}'\n\
634 if @var{op} is @qcode{\"+\"}\n\
635 \n\
636 @item\n\
637 @var{R1}'*@var{R1} = @var{R}'*@var{R} - @var{u}*@var{u}'\n\
638 if @var{op} is @qcode{\"-\"}\n\
639 @end itemize\n\
640 \n\
641 If @var{op} is @qcode{\"-\"}, @var{info} is set to\n\
642 \n\
643 @itemize\n\
644 @item 0 if the downdate was successful,\n\
645 \n\
646 @item 1 if @var{R}'*@var{R} - @var{u}*@var{u}' is not positive definite,\n\
647 \n\
648 @item 2 if @var{R} is singular.\n\
649 @end itemize\n\
650 \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\
653 @end deftypefn")
654 {
655  octave_idx_type nargin = args.length ();
656 
657  octave_value_list retval;
658 
659  if (nargin > 3 || nargin < 2)
660  {
661  print_usage ();
662  return retval;
663  }
664 
665  octave_value argr = args(0);
666  octave_value argu = args(1);
667 
668  if (argr.is_numeric_type () && argu.is_numeric_type ()
669  && (nargin < 3 || args(2).is_string ()))
670  {
671  octave_idx_type n = argr.rows ();
672 
673  std::string op = (nargin < 3) ? "+" : args(2).string_value ();
674 
675  bool down = op == "-";
676 
677  if (down || op == "+")
678  if (argr.columns () == n && argu.rows () == n && argu.columns () == 1)
679  {
680  int err = 0;
681  if (argr.is_single_type () || argu.is_single_type ())
682  {
683  if (argr.is_real_type () && argu.is_real_type ())
684  {
685  // real case
686  FloatMatrix R = argr.float_matrix_value ();
688 
689  FloatCHOL fact;
690  fact.set (R);
691 
692  if (down)
693  err = fact.downdate (u);
694  else
695  fact.update (u);
696 
697  retval(0) = get_chol_r (fact);
698  }
699  else
700  {
701  // complex case
705 
706  FloatComplexCHOL fact;
707  fact.set (R);
708 
709  if (down)
710  err = fact.downdate (u);
711  else
712  fact.update (u);
713 
714  retval(0) = get_chol_r (fact);
715  }
716  }
717  else
718  {
719  if (argr.is_real_type () && argu.is_real_type ())
720  {
721  // real case
722  Matrix R = argr.matrix_value ();
723  ColumnVector u = argu.column_vector_value ();
724 
725  CHOL fact;
726  fact.set (R);
727 
728  if (down)
729  err = fact.downdate (u);
730  else
731  fact.update (u);
732 
733  retval(0) = get_chol_r (fact);
734  }
735  else
736  {
737  // complex case
740 
741  ComplexCHOL fact;
742  fact.set (R);
743 
744  if (down)
745  err = fact.downdate (u);
746  else
747  fact.update (u);
748 
749  retval(0) = get_chol_r (fact);
750  }
751  }
752 
753  if (nargout > 1)
754  retval(1) = err;
755  else if (err == 1)
756  error ("cholupdate: downdate violates positiveness");
757  else if (err == 2)
758  error ("cholupdate: singular matrix");
759  }
760  else
761  error ("cholupdate: dimension mismatch between R and U");
762  else
763  error ("cholupdate: OP must be \"+\" or \"-\"");
764  }
765  else
766  print_usage ();
767 
768  return retval;
769 }
770 
771 /*
772 %!shared A, u, Ac, uc
773 %! A = [ 0.436997 -0.131721 0.124120 -0.061673 ;
774 %! -0.131721 0.738529 0.019851 -0.140295 ;
775 %! 0.124120 0.019851 0.354879 -0.059472 ;
776 %! -0.061673 -0.140295 -0.059472 0.600939 ];
777 %!
778 %! u = [ 0.98950 ;
779 %! 0.39844 ;
780 %! 0.63484 ;
781 %! 0.13351 ];
782 %! Ac = [ 0.5585528 + 0.0000000i -0.1662088 - 0.0315341i 0.0107873 + 0.0236411i -0.0276775 - 0.0186073i ;
783 %! -0.1662088 + 0.0315341i 0.6760061 + 0.0000000i 0.0011452 - 0.0475528i 0.0145967 + 0.0247641i ;
784 %! 0.0107873 - 0.0236411i 0.0011452 + 0.0475528i 0.6263149 - 0.0000000i -0.1585837 - 0.0719763i ;
785 %! -0.0276775 + 0.0186073i 0.0145967 - 0.0247641i -0.1585837 + 0.0719763i 0.6034234 - 0.0000000i ];
786 %!
787 %! uc = [ 0.54267 + 0.91519i ;
788 %! 0.99647 + 0.43141i ;
789 %! 0.83760 + 0.68977i ;
790 %! 0.39160 + 0.90378i ];
791 
792 %!test
793 %! R = chol (A);
794 %! R1 = cholupdate (R, u);
795 %! assert (norm (triu (R1)-R1, Inf), 0);
796 %! assert (norm (R1'*R1 - R'*R - u*u', Inf) < 1e1*eps);
797 %!
798 %! R1 = cholupdate (R1, u, "-");
799 %! assert (norm (triu (R1)-R1, Inf), 0);
800 %! assert (norm (R1 - R, Inf) < 1e1*eps);
801 
802 %!test
803 %! R = chol (Ac);
804 %! R1 = cholupdate (R, uc);
805 %! assert (norm (triu (R1)-R1, Inf), 0);
806 %! assert (norm (R1'*R1 - R'*R - uc*uc', Inf) < 1e1*eps);
807 %!
808 %! R1 = cholupdate (R1, uc, "-");
809 %! assert (norm (triu (R1)-R1, Inf), 0);
810 %! assert (norm (R1 - R, Inf) < 1e1*eps);
811 
812 %!test
813 %! R = chol (single (A));
814 %! R1 = cholupdate (R, single (u));
815 %! assert (norm (triu (R1)-R1, Inf), single (0));
816 %! assert (norm (R1'*R1 - R'*R - single (u*u'), Inf) < 1e1*eps ("single"));
817 %!
818 %! R1 = cholupdate (R1, single (u), "-");
819 %! assert (norm (triu (R1)-R1, Inf), single (0));
820 %! assert (norm (R1 - R, Inf) < 2e1*eps ("single"));
821 
822 %!test
823 %! R = chol (single (Ac));
824 %! R1 = cholupdate (R, single (uc));
825 %! assert (norm (triu (R1)-R1, Inf), single (0));
826 %! assert (norm (R1'*R1 - R'*R - single (uc*uc'), Inf) < 1e1*eps ("single"));
827 %!
828 %! R1 = cholupdate (R1, single (uc), "-");
829 %! assert (norm (triu (R1)-R1, Inf), single (0));
830 %! assert (norm (R1 - R, Inf) < 2e1*eps ("single"));
831 */
832 
833 DEFUN_DLD (cholinsert, args, nargout,
834  "-*- texinfo -*-\n\
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\
843 \n\
844 @itemize\n\
845 @item 0 if the insertion was successful,\n\
846 \n\
847 @item 1 if @var{A1} is not positive definite,\n\
848 \n\
849 @item 2 if @var{R} is singular.\n\
850 @end itemize\n\
851 \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\
854 @end deftypefn")
855 {
856  octave_idx_type nargin = args.length ();
857 
858  octave_value_list retval;
859 
860  if (nargin != 3)
861  {
862  print_usage ();
863  return retval;
864  }
865 
866  octave_value argr = args(0);
867  octave_value argj = args(1);
868  octave_value argu = args(2);
869 
870  if (argr.is_numeric_type () && argu.is_numeric_type ()
871  && argj.is_real_scalar ())
872  {
873  octave_idx_type n = argr.rows ();
874  octave_idx_type j = argj.scalar_value ();
875 
876  if (argr.columns () == n && argu.rows () == n+1 && argu.columns () == 1)
877  {
878  if (j > 0 && j <= n+1)
879  {
880  int err = 0;
881  if (argr.is_single_type () || argu.is_single_type ())
882  {
883  if (argr.is_real_type () && argu.is_real_type ())
884  {
885  // real case
886  FloatMatrix R = argr.float_matrix_value ();
888 
889  FloatCHOL fact;
890  fact.set (R);
891  err = fact.insert_sym (u, j-1);
892 
893  retval(0) = get_chol_r (fact);
894  }
895  else
896  {
897  // complex case
901 
902  FloatComplexCHOL fact;
903  fact.set (R);
904  err = fact.insert_sym (u, j-1);
905 
906  retval(0) = get_chol_r (fact);
907  }
908  }
909  else
910  {
911  if (argr.is_real_type () && argu.is_real_type ())
912  {
913  // real case
914  Matrix R = argr.matrix_value ();
915  ColumnVector u = argu.column_vector_value ();
916 
917  CHOL fact;
918  fact.set (R);
919  err = fact.insert_sym (u, j-1);
920 
921  retval(0) = get_chol_r (fact);
922  }
923  else
924  {
925  // complex case
929 
930  ComplexCHOL fact;
931  fact.set (R);
932  err = fact.insert_sym (u, j-1);
933 
934  retval(0) = get_chol_r (fact);
935  }
936  }
937 
938  if (nargout > 1)
939  retval(1) = err;
940  else if (err == 1)
941  error ("cholinsert: insertion violates positiveness");
942  else if (err == 2)
943  error ("cholinsert: singular matrix");
944  else if (err == 3)
945  error ("cholinsert: diagonal element must be real");
946  }
947  else
948  error ("cholinsert: index J out of range");
949  }
950  else
951  error ("cholinsert: dimension mismatch between R and U");
952  }
953  else
954  print_usage ();
955 
956  return retval;
957 }
958 
959 /*
960 %!test
961 %! u2 = [ 0.35080 ;
962 %! 0.63930 ;
963 %! 3.31057 ;
964 %! -0.13825 ;
965 %! 0.45266 ];
966 %!
967 %! R = chol (A);
968 %!
969 %! j = 3; p = [1:j-1, j+1:5];
970 %! R1 = cholinsert (R, j, u2);
971 %! A1 = R1'*R1;
972 %!
973 %! assert (norm (triu (R1)-R1, Inf), 0);
974 %! assert (norm (A1(p,p) - A, Inf) < 1e1*eps);
975 
976 %!test
977 %! u2 = [ 0.35080 + 0.04298i;
978 %! 0.63930 + 0.23778i;
979 %! 3.31057 + 0.00000i;
980 %! -0.13825 + 0.19879i;
981 %! 0.45266 + 0.50020i];
982 %!
983 %! R = chol (Ac);
984 %!
985 %! j = 3; p = [1:j-1, j+1:5];
986 %! R1 = cholinsert (R, j, u2);
987 %! A1 = R1'*R1;
988 %!
989 %! assert (norm (triu (R1)-R1, Inf), 0);
990 %! assert (norm (A1(p,p) - Ac, Inf) < 1e1*eps);
991 
992 %!test
993 %! u2 = single ([ 0.35080 ;
994 %! 0.63930 ;
995 %! 3.31057 ;
996 %! -0.13825 ;
997 %! 0.45266 ]);
998 %!
999 %! R = chol (single (A));
1000 %!
1001 %! j = 3; p = [1:j-1, j+1:5];
1002 %! R1 = cholinsert (R, j, u2);
1003 %! A1 = R1'*R1;
1004 %!
1005 %! assert (norm (triu (R1)-R1, Inf), single (0));
1006 %! assert (norm (A1(p,p) - A, Inf) < 1e1*eps ("single"));
1007 
1008 %!test
1009 %! u2 = single ([ 0.35080 + 0.04298i;
1010 %! 0.63930 + 0.23778i;
1011 %! 3.31057 + 0.00000i;
1012 %! -0.13825 + 0.19879i;
1013 %! 0.45266 + 0.50020i]);
1014 %!
1015 %! R = chol (single (Ac));
1016 %!
1017 %! j = 3; p = [1:j-1, j+1:5];
1018 %! R1 = cholinsert (R, j, u2);
1019 %! A1 = R1'*R1;
1020 %!
1021 %! assert (norm (triu (R1)-R1, Inf), single (0));
1022 %! assert (norm (A1(p,p) - single (Ac), Inf) < 2e1*eps ("single"));
1023 
1024 %!test
1025 %! cu = chol (triu (A), "upper");
1026 %! cl = chol (tril (A), "lower");
1027 %! assert (cu, cl', eps);
1028 
1029 %!test
1030 %! cca = chol (Ac);
1031 %!
1032 %! ccal = chol (Ac, "lower");
1033 %! ccal2 = chol (tril (Ac), "lower");
1034 %!
1035 %! ccau = chol (Ac, "upper");
1036 %! ccau2 = chol (triu (Ac), "upper");
1037 %!
1038 %! assert (cca'*cca, Ac, eps);
1039 %! assert (ccau'*ccau, Ac, eps);
1040 %! assert (ccau2'*ccau2, Ac, eps);
1041 %!
1042 %! assert (cca, ccal', eps);
1043 %! assert (cca, ccau, eps);
1044 %! assert (cca, ccal2', eps);
1045 %! assert (cca, ccau2, eps);
1046 
1047 %!test
1048 %! cca = chol (single (Ac));
1049 %!
1050 %! ccal = chol (single (Ac), "lower");
1051 %! ccal2 = chol (tril (single (Ac)), "lower");
1052 %!
1053 %! ccau = chol (single (Ac), "upper");
1054 %! ccau2 = chol (triu (single (Ac)), "upper");
1055 %!
1056 %! assert (cca'*cca, single (Ac), eps ("single"));
1057 %! assert (ccau'*ccau, single (Ac), eps ("single"));
1058 %! assert (ccau2'*ccau2, single (Ac), eps ("single"));
1059 %!
1060 %! assert (cca, ccal', eps ("single"));
1061 %! assert (cca, ccau, eps ("single"));
1062 %! assert (cca, ccal2', eps ("single"));
1063 %! assert (cca, ccau2, eps ("single"));
1064 
1065 %!test
1066 %! a = [12, 2, 3, 4;
1067 %! 2, 14, 5, 3;
1068 %! 3, 5, 16, 6;
1069 %! 4, 3, 6, 16];
1070 %!
1071 %! b = [0, 1, 2, 3;
1072 %! -1, 0, 1, 2;
1073 %! -2, -1, 0, 1;
1074 %! -3, -2, -1, 0];
1075 %!
1076 %! ca = a + i*b;
1077 %!
1078 %! cca = chol (ca);
1079 %!
1080 %! ccal = chol (ca, "lower");
1081 %! ccal2 = chol (tril (ca), "lower");
1082 %!
1083 %! ccau = chol (ca, "upper");
1084 %! ccau2 = chol (triu (ca), "upper");
1085 %!
1086 %! assert (cca'*cca, ca, 16*eps);
1087 %! assert (ccau'*ccau, ca, 16*eps);
1088 %! assert (ccau2'*ccau2, ca, 16*eps);
1089 %!
1090 %! assert (cca, ccal', 16*eps);
1091 %! assert (cca, ccau, 16*eps);
1092 %! assert (cca, ccal2', 16*eps);
1093 %! assert (cca, ccau2, 16*eps);
1094 */
1095 
1096 DEFUN_DLD (choldelete, args, ,
1097  "-*- texinfo -*-\n\
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\
1104 @end deftypefn")
1105 {
1106  octave_idx_type nargin = args.length ();
1107 
1108  octave_value_list retval;
1109 
1110  if (nargin != 2)
1111  {
1112  print_usage ();
1113  return retval;
1114  }
1115 
1116  octave_value argr = args(0);
1117  octave_value argj = args(1);
1118 
1119  if (argr.is_numeric_type () && argj.is_real_scalar ())
1120  {
1121  octave_idx_type n = argr.rows ();
1122  octave_idx_type j = argj.scalar_value ();
1123 
1124  if (argr.columns () == n)
1125  {
1126  if (j > 0 && j <= n)
1127  {
1128  if (argr.is_single_type ())
1129  {
1130  if (argr.is_real_type ())
1131  {
1132  // real case
1133  FloatMatrix R = argr.float_matrix_value ();
1134 
1135  FloatCHOL fact;
1136  fact.set (R);
1137  fact.delete_sym (j-1);
1138 
1139  retval(0) = get_chol_r (fact);
1140  }
1141  else
1142  {
1143  // complex case
1145 
1146  FloatComplexCHOL fact;
1147  fact.set (R);
1148  fact.delete_sym (j-1);
1149 
1150  retval(0) = get_chol_r (fact);
1151  }
1152  }
1153  else
1154  {
1155  if (argr.is_real_type ())
1156  {
1157  // real case
1158  Matrix R = argr.matrix_value ();
1159 
1160  CHOL fact;
1161  fact.set (R);
1162  fact.delete_sym (j-1);
1163 
1164  retval(0) = get_chol_r (fact);
1165  }
1166  else
1167  {
1168  // complex case
1169  ComplexMatrix R = argr.complex_matrix_value ();
1170 
1171  ComplexCHOL fact;
1172  fact.set (R);
1173  fact.delete_sym (j-1);
1174 
1175  retval(0) = get_chol_r (fact);
1176  }
1177  }
1178  }
1179  else
1180  error ("choldelete: index J out of range");
1181  }
1182  else
1183  error ("choldelete: matrix R must be square");
1184  }
1185  else
1186  print_usage ();
1187 
1188  return retval;
1189 }
1190 
1191 /*
1192 %!test
1193 %! R = chol (A);
1194 %!
1195 %! j = 3; p = [1:j-1,j+1:4];
1196 %! R1 = choldelete (R, j);
1197 %!
1198 %! assert (norm (triu (R1)-R1, Inf), 0);
1199 %! assert (norm (R1'*R1 - A(p,p), Inf) < 1e1*eps);
1200 
1201 %!test
1202 %! R = chol (Ac);
1203 %!
1204 %! j = 3; p = [1:j-1,j+1:4];
1205 %! R1 = choldelete (R, j);
1206 %!
1207 %! assert (norm (triu (R1)-R1, Inf), 0);
1208 %! assert (norm (R1'*R1 - Ac(p,p), Inf) < 1e1*eps);
1209 
1210 %!test
1211 %! R = chol (single (A));
1212 %!
1213 %! j = 3; p = [1:j-1,j+1:4];
1214 %! R1 = choldelete (R, j);
1215 %!
1216 %! assert (norm (triu (R1)-R1, Inf), single (0));
1217 %! assert (norm (R1'*R1 - single (A(p,p)), Inf) < 1e1*eps ("single"));
1218 
1219 %!test
1220 %! R = chol (single (Ac));
1221 %!
1222 %! j = 3; p = [1:j-1,j+1:4];
1223 %! R1 = choldelete (R,j);
1224 %!
1225 %! assert (norm (triu (R1)-R1, Inf), single (0));
1226 %! assert (norm (R1'*R1 - single (Ac(p,p)), Inf) < 1e1*eps ("single"));
1227 */
1228 
1229 DEFUN_DLD (cholshift, args, ,
1230  "-*- texinfo -*-\n\
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\
1237  or @*\n\
1238 @code{p = [1:j-1, shift(j:i,-1), i+1:n]} if @w{@var{j} < @var{i}}. @*\n\
1239 \n\
1240 @seealso{chol, cholupdate, cholinsert, choldelete}\n\
1241 @end deftypefn")
1242 {
1243  octave_idx_type nargin = args.length ();
1244 
1245  octave_value_list retval;
1246 
1247  if (nargin != 3)
1248  {
1249  print_usage ();
1250  return retval;
1251  }
1252 
1253  octave_value argr = args(0);
1254  octave_value argi = args(1);
1255  octave_value argj = args(2);
1256 
1257  if (argr.is_numeric_type ()
1258  && argi.is_real_scalar () && argj.is_real_scalar ())
1259  {
1260  octave_idx_type n = argr.rows ();
1261  octave_idx_type i = argi.scalar_value ();
1262  octave_idx_type j = argj.scalar_value ();
1263 
1264  if (argr.columns () == n)
1265  {
1266  if (j > 0 && j <= n+1 && i > 0 && i <= n+1)
1267  {
1268 
1269  if (argr.is_single_type () && argi.is_single_type () &&
1270  argj.is_single_type ())
1271  {
1272  if (argr.is_real_type ())
1273  {
1274  // real case
1275  FloatMatrix R = argr.float_matrix_value ();
1276 
1277  FloatCHOL fact;
1278  fact.set (R);
1279  fact.shift_sym (i-1, j-1);
1280 
1281  retval(0) = get_chol_r (fact);
1282  }
1283  else
1284  {
1285  // complex case
1287 
1288  FloatComplexCHOL fact;
1289  fact.set (R);
1290  fact.shift_sym (i-1, j-1);
1291 
1292  retval(0) = get_chol_r (fact);
1293  }
1294  }
1295  else
1296  {
1297  if (argr.is_real_type ())
1298  {
1299  // real case
1300  Matrix R = argr.matrix_value ();
1301 
1302  CHOL fact;
1303  fact.set (R);
1304  fact.shift_sym (i-1, j-1);
1305 
1306  retval(0) = get_chol_r (fact);
1307  }
1308  else
1309  {
1310  // complex case
1311  ComplexMatrix R = argr.complex_matrix_value ();
1312 
1313  ComplexCHOL fact;
1314  fact.set (R);
1315  fact.shift_sym (i-1, j-1);
1316 
1317  retval(0) = get_chol_r (fact);
1318  }
1319  }
1320  }
1321  else
1322  error ("cholshift: index I or J is out of range");
1323  }
1324  else
1325  error ("cholshift: R must be a square matrix");
1326  }
1327  else
1328  print_usage ();
1329 
1330  return retval;
1331 }
1332 
1333 /*
1334 %!test
1335 %! R = chol (A);
1336 %!
1337 %! i = 1; j = 3; p = [1:i-1, shift(i:j,-1), j+1:4];
1338 %! R1 = cholshift (R, i, j);
1339 %!
1340 %! assert (norm (triu (R1)-R1, Inf), 0);
1341 %! assert (norm (R1'*R1 - A(p,p), Inf) < 1e1*eps);
1342 %!
1343 %! j = 1; i = 3; p = [1:j-1, shift(j:i,+1), i+1:4];
1344 %! R1 = cholshift (R, i, j);
1345 %!
1346 %! assert (norm (triu (R1) - R1, Inf), 0);
1347 %! assert (norm (R1'*R1 - A(p,p), Inf) < 1e1*eps);
1348 
1349 %!test
1350 %! R = chol (Ac);
1351 %!
1352 %! i = 1; j = 3; p = [1:i-1, shift(i:j,-1), j+1:4];
1353 %! R1 = cholshift (R, i, j);
1354 %!
1355 %! assert (norm (triu (R1)-R1, Inf), 0);
1356 %! assert (norm (R1'*R1 - Ac(p,p), Inf) < 1e1*eps);
1357 %!
1358 %! j = 1; i = 3; p = [1:j-1, shift(j:i,+1), i+1:4];
1359 %! R1 = cholshift (R, i, j);
1360 %!
1361 %! assert (norm (triu (R1)-R1, Inf), 0);
1362 %! assert (norm (R1'*R1 - Ac(p,p), Inf) < 1e1*eps);
1363 
1364 %!test
1365 %! R = chol (single (A));
1366 %!
1367 %! i = 1; j = 3; p = [1:i-1, shift(i:j,-1), j+1:4];
1368 %! R1 = cholshift (R, i, j);
1369 %!
1370 %! assert (norm (triu (R1)-R1, Inf), 0);
1371 %! assert (norm (R1'*R1 - single (A(p,p)), Inf) < 1e1*eps ("single"));
1372 %!
1373 %! j = 1; i = 3; p = [1:j-1, shift(j:i,+1), i+1:4];
1374 %! R1 = cholshift (R, i, j);
1375 %!
1376 %! assert (norm (triu (R1)-R1, Inf), 0);
1377 %! assert (norm (R1'*R1 - single (A(p,p)), Inf) < 1e1*eps ("single"));
1378 
1379 %!test
1380 %! R = chol (single (Ac));
1381 %!
1382 %! i = 1; j = 3; p = [1:i-1, shift(i:j,-1), j+1:4];
1383 %! R1 = cholshift (R, i, j);
1384 %!
1385 %! assert (norm (triu (R1)-R1, Inf), 0);
1386 %! assert (norm (R1'*R1 - single (Ac(p,p)), Inf) < 1e1*eps ("single"));
1387 %!
1388 %! j = 1; i = 3; p = [1:j-1, shift(j:i,+1), i+1:4];
1389 %! R1 = cholshift (R, i, j);
1390 %!
1391 %! assert (norm (triu (R1)-R1, Inf), 0);
1392 %! assert (norm (R1'*R1 - single (Ac(p,p)), Inf) < 1e1*eps ("single"));
1393 */