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 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178
|
<!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: Pivoted Cholesky Decomposition</title>
<meta name="description" content="GNU Scientific Library – Reference Manual: Pivoted Cholesky Decomposition">
<meta name="keywords" content="GNU Scientific Library – Reference Manual: Pivoted 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="Modified-Cholesky-Decomposition.html#Modified-Cholesky-Decomposition" rel="next" title="Modified Cholesky Decomposition">
<link href="Cholesky-Decomposition.html#Cholesky-Decomposition" rel="previous" title="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="Pivoted-Cholesky-Decomposition"></a>
<div class="header">
<p>
Next: <a href="Modified-Cholesky-Decomposition.html#Modified-Cholesky-Decomposition" accesskey="n" rel="next">Modified Cholesky Decomposition</a>, Previous: <a href="Cholesky-Decomposition.html#Cholesky-Decomposition" accesskey="p" rel="previous">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="Pivoted-Cholesky-Decomposition-1"></a>
<h3 class="section">14.7 Pivoted Cholesky Decomposition</h3>
<a name="index-Cholesky-decomposition_002c-pivoted"></a>
<a name="index-Pivoted-Cholesky-Decomposition"></a>
<p>A symmetric, positive definite square matrix <em>A</em> has an alternate
Cholesky decomposition into a product of a lower unit triangular matrix <em>L</em>,
a diagonal matrix <em>D</em> and <em>L^T</em>, given by <em>L D L^T</em>. This is
equivalent to the Cholesky formulation discussed above, with
the standard Cholesky lower triangular factor given by <em>L D^{1 \over 2}</em>.
For ill-conditioned matrices, it can help to use a pivoting strategy to
prevent the entries of <em>D</em> and <em>L</em> from growing too large, and also
ensure <em>D_1 \ge D_2 \ge \cdots \ge D_n > 0</em>, where <em>D_i</em> are
the diagonal entries of <em>D</em>. The final decomposition is given by
</p>
<div class="example">
<pre class="example">P A P^T = L D L^T
</pre></div>
<p>where <em>P</em> is a permutation matrix.
</p>
<dl>
<dt><a name="index-gsl_005flinalg_005fpcholesky_005fdecomp"></a>Function: <em>int</em> <strong>gsl_linalg_pcholesky_decomp</strong> <em>(gsl_matrix * <var>A</var>, gsl_permutation * <var>p</var>)</em></dt>
<dd><p>This function factors the symmetric, positive-definite square matrix
<var>A</var> into the Pivoted Cholesky decomposition <em>P A 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.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005flinalg_005fpcholesky_005fsolve"></a>Function: <em>int</em> <strong>gsl_linalg_pcholesky_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 system <em>A x = b</em> using the Pivoted Cholesky
decomposition of <em>A</em> held in the matrix <var>LDLT</var> and permutation
<var>p</var> which must have been previously computed by <code>gsl_linalg_pcholesky_decomp</code>.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005flinalg_005fpcholesky_005fsvx"></a>Function: <em>int</em> <strong>gsl_linalg_pcholesky_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 system <em>A x = b</em> in-place using the Pivoted Cholesky
decomposition of <em>A</em> held in the matrix <var>LDLT</var> and permutation
<var>p</var> which must have been previously computed by <code>gsl_linalg_pcholesky_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_005fpcholesky_005fdecomp2"></a>Function: <em>int</em> <strong>gsl_linalg_pcholesky_decomp2</strong> <em>(gsl_matrix * <var>A</var>, gsl_permutation * <var>p</var>, gsl_vector * <var>S</var>)</em></dt>
<dd><p>This function computes the pivoted Cholesky factorization of the matrix
<em>S A S</em>, where the input matrix <var>A</var> is symmetric and positive
definite, and the diagonal scaling matrix <var>S</var> is computed to reduce the
condition number of <var>A</var> as much as possible. See
<a href="Cholesky-Decomposition.html#Cholesky-Decomposition">Cholesky Decomposition</a> for more information on the matrix <var>S</var>.
The Pivoted Cholesky decomposition satisfies <em>P S A S 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 scaling transformation is stored in <var>S</var> on output.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005flinalg_005fpcholesky_005fsolve2"></a>Function: <em>int</em> <strong>gsl_linalg_pcholesky_solve2</strong> <em>(const gsl_matrix * <var>LDLT</var>, const gsl_permutation * <var>p</var>, const gsl_vector * <var>S</var>, const gsl_vector * <var>b</var>, gsl_vector * <var>x</var>)</em></dt>
<dd><p>This function solves the system <em>(S A S) (S^{-1} x) = S b</em> using the Pivoted Cholesky
decomposition of <em>S A S</em> held in the matrix <var>LDLT</var>, permutation
<var>p</var>, and vector <var>S</var>, which must have been previously computed by
<code>gsl_linalg_pcholesky_decomp2</code>.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005flinalg_005fpcholesky_005fsvx2"></a>Function: <em>int</em> <strong>gsl_linalg_pcholesky_svx2</strong> <em>(const gsl_matrix * <var>LDLT</var>, const gsl_permutation * <var>p</var>, const gsl_vector * <var>S</var>, gsl_vector * <var>x</var>)</em></dt>
<dd><p>This function solves the system <em>(S A S) (S^{-1} x) = S b</em> in-place using the Pivoted Cholesky
decomposition of <em>S A S</em> held in the matrix <var>LDLT</var>, permutation
<var>p</var> and vector <var>S</var>, which must have been previously computed by
<code>gsl_linalg_pcholesky_decomp2</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_005fpcholesky_005finvert"></a>Function: <em>int</em> <strong>gsl_linalg_pcholesky_invert</strong> <em>(const gsl_matrix * <var>LDLT</var>, const gsl_permutation * <var>p</var>, gsl_matrix * <var>Ainv</var>)</em></dt>
<dd><p>This function computes the inverse of the matrix <em>A</em>, using the Pivoted
Cholesky decomposition stored in <var>LDLT</var> and <var>p</var>. On output, the
matrix <var>Ainv</var> contains <em>A^{-1}</em>.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005flinalg_005fpcholesky_005frcond"></a>Function: <em>int</em> <strong>gsl_linalg_pcholesky_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 symmetric positive
definite matrix <em>A</em>, using its pivoted Cholesky decomposition provided in <var>LDLT</var>.
The reciprocal condition number estimate, defined as <em>1 / (||A||_1 \cdot ||A^{-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="Modified-Cholesky-Decomposition.html#Modified-Cholesky-Decomposition" accesskey="n" rel="next">Modified Cholesky Decomposition</a>, Previous: <a href="Cholesky-Decomposition.html#Cholesky-Decomposition" accesskey="p" rel="previous">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>
|