1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146
|
<!DOCTYPE html PUBLIC "-//W3C//DTD HTML 4.01 Transitional//EN" "http://www.w3.org/TR/html4/loose.dtd">
<html>
<!-- Copyright (C) 1996, 1997, 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009, 2010, 2011, 2012, 2013, 2014, 2015, 2016 The GSL Team.
Permission is granted to copy, distribute and/or modify this document
under the terms of the GNU Free Documentation License, Version 1.3 or
any later version published by the Free Software Foundation; with the
Invariant Sections being "GNU General Public License" and "Free Software
Needs Free Documentation", the Front-Cover text being "A GNU Manual",
and with the Back-Cover Text being (a) (see below). A copy of the
license is included in the section entitled "GNU Free Documentation
License".
(a) The Back-Cover Text is: "You have the freedom to copy and modify this
GNU Manual." -->
<!-- Created by GNU Texinfo 5.1, http://www.gnu.org/software/texinfo/ -->
<head>
<title>GNU Scientific Library – Reference Manual: Modified Cholesky Decomposition</title>
<meta name="description" content="GNU Scientific Library – Reference Manual: Modified Cholesky Decomposition">
<meta name="keywords" content="GNU Scientific Library – Reference Manual: Modified Cholesky Decomposition">
<meta name="resource-type" content="document">
<meta name="distribution" content="global">
<meta name="Generator" content="makeinfo">
<meta http-equiv="Content-Type" content="text/html; charset=utf-8">
<link href="index.html#Top" rel="start" title="Top">
<link href="Function-Index.html#Function-Index" rel="index" title="Function Index">
<link href="Linear-Algebra.html#Linear-Algebra" rel="up" title="Linear Algebra">
<link href="Tridiagonal-Decomposition-of-Real-Symmetric-Matrices.html#Tridiagonal-Decomposition-of-Real-Symmetric-Matrices" rel="next" title="Tridiagonal Decomposition of Real Symmetric Matrices">
<link href="Pivoted-Cholesky-Decomposition.html#Pivoted-Cholesky-Decomposition" rel="previous" title="Pivoted Cholesky Decomposition">
<style type="text/css">
<!--
a.summary-letter {text-decoration: none}
blockquote.smallquotation {font-size: smaller}
div.display {margin-left: 3.2em}
div.example {margin-left: 3.2em}
div.indentedblock {margin-left: 3.2em}
div.lisp {margin-left: 3.2em}
div.smalldisplay {margin-left: 3.2em}
div.smallexample {margin-left: 3.2em}
div.smallindentedblock {margin-left: 3.2em; font-size: smaller}
div.smalllisp {margin-left: 3.2em}
kbd {font-style:oblique}
pre.display {font-family: inherit}
pre.format {font-family: inherit}
pre.menu-comment {font-family: serif}
pre.menu-preformatted {font-family: serif}
pre.smalldisplay {font-family: inherit; font-size: smaller}
pre.smallexample {font-size: smaller}
pre.smallformat {font-family: inherit; font-size: smaller}
pre.smalllisp {font-size: smaller}
span.nocodebreak {white-space:nowrap}
span.nolinebreak {white-space:nowrap}
span.roman {font-family:serif; font-weight:normal}
span.sansserif {font-family:sans-serif; font-weight:normal}
ul.no-bullet {list-style: none}
-->
</style>
</head>
<body lang="en" bgcolor="#FFFFFF" text="#000000" link="#0000FF" vlink="#800080" alink="#FF0000">
<a name="Modified-Cholesky-Decomposition"></a>
<div class="header">
<p>
Next: <a href="Tridiagonal-Decomposition-of-Real-Symmetric-Matrices.html#Tridiagonal-Decomposition-of-Real-Symmetric-Matrices" accesskey="n" rel="next">Tridiagonal Decomposition of Real Symmetric Matrices</a>, Previous: <a href="Pivoted-Cholesky-Decomposition.html#Pivoted-Cholesky-Decomposition" accesskey="p" rel="previous">Pivoted Cholesky Decomposition</a>, Up: <a href="Linear-Algebra.html#Linear-Algebra" accesskey="u" rel="up">Linear Algebra</a> [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
<hr>
<a name="Modified-Cholesky-Decomposition-1"></a>
<h3 class="section">14.8 Modified Cholesky Decomposition</h3>
<a name="index-Cholesky-decomposition_002c-modified"></a>
<a name="index-Modified-Cholesky-Decomposition"></a>
<p>The modified Cholesky decomposition is suitable for solving systems
<em>A x = b</em> where <em>A</em> is a symmetric indefinite matrix. Such
matrices arise in nonlinear optimization algorithms. The standard
Cholesky decomposition requires a positive definite matrix and would
fail in this case. Instead of resorting to a method like QR or SVD,
which do not take into account the symmetry of the matrix, we can
instead introduce a small perturbation to the matrix <em>A</em> to
make it positive definite, and then use a Cholesky decomposition on
the perturbed matrix. The resulting decomposition satisfies
</p>
<div class="example">
<pre class="example">P (A + E) P^T = L D L^T
</pre></div>
<p>where <em>P</em> is a permutation matrix, <em>E</em> is a diagonal
perturbation matrix, <em>L</em> is unit lower triangular, and
<em>D</em> is diagonal. If <em>A</em> is sufficiently positive
definite, then the perturbation matrix <em>E</em> will be zero
and this method is equivalent to the pivoted Cholesky algorithm.
For indefinite matrices, the perturbation matrix <em>E</em> is
computed to ensure that <em>A + E</em> is positive definite and
well conditioned.
</p>
<dl>
<dt><a name="index-gsl_005flinalg_005fmcholesky_005fdecomp"></a>Function: <em>int</em> <strong>gsl_linalg_mcholesky_decomp</strong> <em>(gsl_matrix * <var>A</var>, gsl_permutation * <var>p</var>, gsl_vector * <var>E</var>)</em></dt>
<dd><p>This function factors the symmetric, indefinite square matrix
<var>A</var> into the Modified Cholesky decomposition <em>P (A + E) P^T = L D L^T</em>.
On input, the values from the diagonal and lower-triangular part of the matrix <var>A</var> are
used to construct the factorization. On output the diagonal of the input matrix <var>A</var> stores the diagonal
elements of <em>D</em>, and the lower triangular portion of <var>A</var>
contains the matrix <em>L</em>. Since <em>L</em> has ones on its diagonal these
do not need to be explicitely stored. The upper triangular portion of <var>A</var>
is unmodified. The permutation matrix <em>P</em> is
stored in <var>p</var> on output. The diagonal perturbation matrix is stored in
<var>E</var> on output. The parameter <var>E</var> may be set to NULL if it is not
required.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005flinalg_005fmcholesky_005fsolve"></a>Function: <em>int</em> <strong>gsl_linalg_mcholesky_solve</strong> <em>(const gsl_matrix * <var>LDLT</var>, const gsl_permutation * <var>p</var>, const gsl_vector * <var>b</var>, gsl_vector * <var>x</var>)</em></dt>
<dd><p>This function solves the perturbed system <em>(A + E) x = b</em> using the Cholesky
decomposition of <em>A + E</em> held in the matrix <var>LDLT</var> and permutation
<var>p</var> which must have been previously computed by <code>gsl_linalg_mcholesky_decomp</code>.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005flinalg_005fmcholesky_005fsvx"></a>Function: <em>int</em> <strong>gsl_linalg_mcholesky_svx</strong> <em>(const gsl_matrix * <var>LDLT</var>, const gsl_permutation * <var>p</var>, gsl_vector * <var>x</var>)</em></dt>
<dd><p>This function solves the perturbed system <em>(A + E) x = b</em> in-place using the Cholesky
decomposition of <em>A + E</em> held in the matrix <var>LDLT</var> and permutation
<var>p</var> which must have been previously computed by <code>gsl_linalg_mcholesky_decomp</code>.
On input, <var>x</var> contains the right hand side vector <em>b</em> which is
replaced by the solution vector on output.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005flinalg_005fmcholesky_005frcond"></a>Function: <em>int</em> <strong>gsl_linalg_mcholesky_rcond</strong> <em>(const gsl_matrix * <var>LDLT</var>, const gsl_permutation * <var>p</var>, double * <var>rcond</var>, gsl_vector * <var>work</var>)</em></dt>
<dd><p>This function estimates the reciprocal condition number (using the 1-norm) of the perturbed matrix
<em>A + E</em>, using its pivoted Cholesky decomposition provided in <var>LDLT</var>.
The reciprocal condition number estimate, defined as <em>1 / (||A + E||_1 \cdot ||(A + E)^{-1}||_1)</em>, is stored
in <var>rcond</var>. Additional workspace of size <em>3 N</em> is required in <var>work</var>.
</p></dd></dl>
<hr>
<div class="header">
<p>
Next: <a href="Tridiagonal-Decomposition-of-Real-Symmetric-Matrices.html#Tridiagonal-Decomposition-of-Real-Symmetric-Matrices" accesskey="n" rel="next">Tridiagonal Decomposition of Real Symmetric Matrices</a>, Previous: <a href="Pivoted-Cholesky-Decomposition.html#Pivoted-Cholesky-Decomposition" accesskey="p" rel="previous">Pivoted Cholesky Decomposition</a>, Up: <a href="Linear-Algebra.html#Linear-Algebra" accesskey="u" rel="up">Linear Algebra</a> [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
</body>
</html>
|