File: Cholesky-Decomposition.html

package info (click to toggle)
gsl-ref-html 2.3-1
  • links: PTS
  • area: non-free
  • in suites: bullseye, buster, sid
  • size: 6,876 kB
  • ctags: 4,574
  • sloc: makefile: 35
file content (227 lines) | stat: -rw-r--r-- 13,479 bytes parent folder | download
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
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
<!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 &ndash; Reference Manual: Cholesky Decomposition</title>

<meta name="description" content="GNU Scientific Library &ndash; Reference Manual: Cholesky Decomposition">
<meta name="keywords" content="GNU Scientific Library &ndash; Reference Manual: 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="Pivoted-Cholesky-Decomposition.html#Pivoted-Cholesky-Decomposition" rel="next" title="Pivoted Cholesky Decomposition">
<link href="Singular-Value-Decomposition.html#Singular-Value-Decomposition" rel="previous" title="Singular Value 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="Cholesky-Decomposition"></a>
<div class="header">
<p>
Next: <a href="Pivoted-Cholesky-Decomposition.html#Pivoted-Cholesky-Decomposition" accesskey="n" rel="next">Pivoted Cholesky Decomposition</a>, Previous: <a href="Singular-Value-Decomposition.html#Singular-Value-Decomposition" accesskey="p" rel="previous">Singular Value Decomposition</a>, Up: <a href="Linear-Algebra.html#Linear-Algebra" accesskey="u" rel="up">Linear Algebra</a> &nbsp; [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
<hr>
<a name="Cholesky-Decomposition-1"></a>
<h3 class="section">14.6 Cholesky Decomposition</h3>
<a name="index-Cholesky-decomposition"></a>
<a name="index-square-root-of-a-matrix_002c-Cholesky-decomposition"></a>
<a name="index-matrix-square-root_002c-Cholesky-decomposition"></a>

<p>A symmetric, positive definite square matrix <em>A</em> has a Cholesky
decomposition into a product of a lower triangular matrix <em>L</em> and
its transpose <em>L^T</em>,
</p>
<div class="example">
<pre class="example">A = L L^T
</pre></div>

<p>This is sometimes referred to as taking the square-root of a matrix. The
Cholesky decomposition can only be carried out when all the eigenvalues
of the matrix are positive.  This decomposition can be used to convert
the linear system <em>A x = b</em> into a pair of triangular systems
(<em>L y = b</em>, <em>L^T x = y</em>), which can be solved by forward and
back-substitution.
</p>
<p>If the matrix <em>A</em> is near singular, it is sometimes possible to reduce
the condition number and recover a more accurate solution vector <em>x</em>
by scaling as
</p>
<div class="example">
<pre class="example">( S A S ) ( S^(-1) x ) = S b
</pre></div>

<p>where <em>S</em> is a diagonal matrix whose elements are given by
<em>S_{ii} = 1/\sqrt{A_{ii}}</em>. This scaling is also known as
Jacobi preconditioning. There are routines below to solve
both the scaled and unscaled systems.
</p>
<dl>
<dt><a name="index-gsl_005flinalg_005fcholesky_005fdecomp1"></a>Function: <em>int</em> <strong>gsl_linalg_cholesky_decomp1</strong> <em>(gsl_matrix * <var>A</var>)</em></dt>
<dt><a name="index-gsl_005flinalg_005fcomplex_005fcholesky_005fdecomp"></a>Function: <em>int</em> <strong>gsl_linalg_complex_cholesky_decomp</strong> <em>(gsl_matrix_complex * <var>A</var>)</em></dt>
<dd><p>These functions factorize the symmetric, positive-definite square matrix
<var>A</var> into the Cholesky decomposition <em>A = L L^T</em> (or
<em>A = L L^H</em>
for the complex case). On input, the values from the diagonal and lower-triangular
part of the matrix <var>A</var> are used (the upper triangular part is ignored).  On output
the diagonal and lower triangular part of the input matrix <var>A</var> contain the matrix
<em>L</em>, while the upper triangular part is unmodified.  If the matrix is not
positive-definite then the decomposition will fail, returning the
error code <code>GSL_EDOM</code>.  
</p>
<p>When testing whether a matrix is positive-definite, disable the error
handler first to avoid triggering an error.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fcholesky_005fdecomp"></a>Function: <em>int</em> <strong>gsl_linalg_cholesky_decomp</strong> <em>(gsl_matrix * <var>A</var>)</em></dt>
<dd><p>This function is now deprecated and is provided only for backward compatibility.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fcholesky_005fsolve"></a>Function: <em>int</em> <strong>gsl_linalg_cholesky_solve</strong> <em>(const gsl_matrix * <var>cholesky</var>, const gsl_vector * <var>b</var>, gsl_vector * <var>x</var>)</em></dt>
<dt><a name="index-gsl_005flinalg_005fcomplex_005fcholesky_005fsolve"></a>Function: <em>int</em> <strong>gsl_linalg_complex_cholesky_solve</strong> <em>(const gsl_matrix_complex * <var>cholesky</var>, const gsl_vector_complex * <var>b</var>, gsl_vector_complex * <var>x</var>)</em></dt>
<dd><p>These functions solve the system <em>A x = b</em> using the Cholesky
decomposition of <em>A</em> held in the matrix <var>cholesky</var> which must
have been previously computed by <code>gsl_linalg_cholesky_decomp</code> or
<code>gsl_linalg_complex_cholesky_decomp</code>.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fcholesky_005fsvx"></a>Function: <em>int</em> <strong>gsl_linalg_cholesky_svx</strong> <em>(const gsl_matrix * <var>cholesky</var>, gsl_vector * <var>x</var>)</em></dt>
<dt><a name="index-gsl_005flinalg_005fcomplex_005fcholesky_005fsvx"></a>Function: <em>int</em> <strong>gsl_linalg_complex_cholesky_svx</strong> <em>(const gsl_matrix_complex * <var>cholesky</var>, gsl_vector_complex * <var>x</var>)</em></dt>
<dd><p>These functions solve the system <em>A x = b</em> in-place using the
Cholesky decomposition of <em>A</em> held in the matrix <var>cholesky</var>
which must have been previously computed by
<code>gsl_linalg_cholesky_decomp</code> or
<code>gsl_linalg_complex_cholesky_decomp</code>. On input <var>x</var> should
contain the right-hand side <em>b</em>, which is replaced by the
solution on output.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fcholesky_005finvert"></a>Function: <em>int</em> <strong>gsl_linalg_cholesky_invert</strong> <em>(gsl_matrix * <var>cholesky</var>)</em></dt>
<dt><a name="index-gsl_005flinalg_005fcomplex_005fcholesky_005finvert"></a>Function: <em>int</em> <strong>gsl_linalg_complex_cholesky_invert</strong> <em>(gsl_matrix_complex * <var>cholesky</var>)</em></dt>
<dd><p>These functions compute the inverse of a matrix from its Cholesky
decomposition <var>cholesky</var>, which must have been previously computed
by <code>gsl_linalg_cholesky_decomp</code> or
<code>gsl_linalg_complex_cholesky_decomp</code>.  On output, the inverse is
stored in-place in <var>cholesky</var>.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fcholesky_005fdecomp2"></a>Function: <em>int</em> <strong>gsl_linalg_cholesky_decomp2</strong> <em>(gsl_matrix * <var>A</var>, gsl_vector * <var>S</var>)</em></dt>
<dd><p>This function calculates a diagonal scaling transformation <em>S</em> for
the symmetric, positive-definite square matrix <var>A</var>, and then
computes the Cholesky decomposition <em>S A S = L L^T</em>.
On input, the values from the diagonal and lower-triangular part of the matrix <var>A</var> are
used (the upper triangular part is ignored).  On output the diagonal and lower triangular part
of the input matrix <var>A</var> contain the matrix <em>L</em>, while the upper triangular part
of the input matrix is overwritten with <em>L^T</em> (the diagonal terms being
identical for both <em>L</em> and <em>L^T</em>).  If the matrix is not
positive-definite then the decomposition will fail, returning the
error code <code>GSL_EDOM</code>. The diagonal scale factors are stored in <var>S</var>
on output.
</p>
<p>When testing whether a matrix is positive-definite, disable the error
handler first to avoid triggering an error.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fcholesky_005fsolve2"></a>Function: <em>int</em> <strong>gsl_linalg_cholesky_solve2</strong> <em>(const gsl_matrix * <var>cholesky</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 Cholesky
decomposition of <em>S A S</em> held in the matrix <var>cholesky</var> which must
have been previously computed by <code>gsl_linalg_cholesky_decomp2</code>.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fcholesky_005fsvx2"></a>Function: <em>int</em> <strong>gsl_linalg_cholesky_svx2</strong> <em>(const gsl_matrix * <var>cholesky</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
Cholesky decomposition of <em>S A S</em> held in the matrix <var>cholesky</var>
which must have been previously computed by
<code>gsl_linalg_cholesky_decomp2</code>.  On input <var>x</var> should
contain the right-hand side <em>b</em>, which is replaced by the
solution on output.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fcholesky_005fscale"></a>Function: <em>int</em> <strong>gsl_linalg_cholesky_scale</strong> <em>(const gsl_matrix * <var>A</var>, gsl_vector * <var>S</var>)</em></dt>
<dd><p>This function calculates a diagonal scaling transformation of the
symmetric, positive definite matrix <var>A</var>, such that
<em>S A S</em> has a condition number within a factor of <em>N</em>
of the matrix of smallest possible condition number over all
possible diagonal scalings. On output, <var>S</var> contains the
scale factors, given by <em>S_i = 1/\sqrt{A_{ii}}</em>.
For any <em>A_{ii} \le 0</em>, the corresponding scale factor <em>S_i</em>
is set to <em>1</em>.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fcholesky_005fscale_005fapply"></a>Function: <em>int</em> <strong>gsl_linalg_cholesky_scale_apply</strong> <em>(gsl_matrix * <var>A</var>, const gsl_vector * <var>S</var>)</em></dt>
<dd><p>This function applies the scaling transformation <var>S</var> to the matrix <var>A</var>. On output,
<var>A</var> is replaced by <em>S A S</em>.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fcholesky_005frcond"></a>Function: <em>int</em> <strong>gsl_linalg_cholesky_rcond</strong> <em>(const gsl_matrix * <var>cholesky</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 Cholesky decomposition provided in <var>cholesky</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="Pivoted-Cholesky-Decomposition.html#Pivoted-Cholesky-Decomposition" accesskey="n" rel="next">Pivoted Cholesky Decomposition</a>, Previous: <a href="Singular-Value-Decomposition.html#Singular-Value-Decomposition" accesskey="p" rel="previous">Singular Value Decomposition</a>, Up: <a href="Linear-Algebra.html#Linear-Algebra" accesskey="u" rel="up">Linear Algebra</a> &nbsp; [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>



</body>
</html>