File: Cholesky-Decomposition.html

package info (click to toggle)
gsl-ref-html 1.10-1
  • links: PTS
  • area: main
  • in suites: lenny
  • size: 4,496 kB
  • ctags: 2,955
  • sloc: makefile: 33
file content (106 lines) | stat: -rw-r--r-- 6,159 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
<html lang="en">
<head>
<title>Cholesky Decomposition - GNU Scientific Library -- Reference Manual</title>
<meta http-equiv="Content-Type" content="text/html">
<meta name="description" content="GNU Scientific Library -- Reference Manual">
<meta name="generator" content="makeinfo 4.8">
<link title="Top" rel="start" href="index.html#Top">
<link rel="up" href="Linear-Algebra.html" title="Linear Algebra">
<link rel="prev" href="Singular-Value-Decomposition.html" title="Singular Value Decomposition">
<link rel="next" href="Tridiagonal-Decomposition-of-Real-Symmetric-Matrices.html" title="Tridiagonal Decomposition of Real Symmetric Matrices">
<link href="http://www.gnu.org/software/texinfo/" rel="generator-home" title="Texinfo Homepage">
<!--
Copyright (C) 1996, 1997, 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007 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.2 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 freedom to copy and modify this
GNU Manual, like GNU software.''-->
<meta http-equiv="Content-Style-Type" content="text/css">
<style type="text/css"><!--
  pre.display { font-family:inherit }
  pre.format  { font-family:inherit }
  pre.smalldisplay { font-family:inherit; font-size:smaller }
  pre.smallformat  { font-family:inherit; font-size:smaller }
  pre.smallexample { font-size:smaller }
  pre.smalllisp    { font-size:smaller }
  span.sc    { font-variant:small-caps }
  span.roman { font-family:serif; font-weight:normal; } 
  span.sansserif { font-family:sans-serif; font-weight:normal; } 
--></style>
</head>
<body>
<div class="node">
<p>
<a name="Cholesky-Decomposition"></a>
Next:&nbsp;<a rel="next" accesskey="n" href="Tridiagonal-Decomposition-of-Real-Symmetric-Matrices.html">Tridiagonal Decomposition of Real Symmetric Matrices</a>,
Previous:&nbsp;<a rel="previous" accesskey="p" href="Singular-Value-Decomposition.html">Singular Value Decomposition</a>,
Up:&nbsp;<a rel="up" accesskey="u" href="Linear-Algebra.html">Linear Algebra</a>
<hr>
</div>

<h3 class="section">13.5 Cholesky Decomposition</h3>

<p><a name="index-Cholesky-decomposition-1279"></a><a name="index-square-root-of-a-matrix_002c-Cholesky-decomposition-1280"></a><a name="index-matrix-square-root_002c-Cholesky-decomposition-1281"></a>
A symmetric, positive definite square matrix A has a Cholesky
decomposition into a product of a lower triangular matrix L and
its transpose L^T,

<pre class="example">     A = L L^T
</pre>
   <p class="noindent">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 A x = b into a pair of triangular systems
(L y = b, L^T x = y), which can be solved by forward and
back-substitution.

<div class="defun">
&mdash; Function: int <b>gsl_linalg_cholesky_decomp</b> (<var>gsl_matrix * A</var>)<var><a name="index-gsl_005flinalg_005fcholesky_005fdecomp-1282"></a></var><br>
&mdash; Function: int <b>gsl_linalg_complex_cholesky_decomp</b> (<var>gsl_matrix_complex * A</var>)<var><a name="index-gsl_005flinalg_005fcomplex_005fcholesky_005fdecomp-1283"></a></var><br>
<blockquote><p>These functions factorize the positive-definite square matrix
<var>A</var> into the Cholesky decomposition A = L L^T (or
<!-- {$A = L L^{\dagger}$} -->
A = L L^H
for the complex case). On input
only the diagonal and lower-triangular part of the matrix <var>A</var> are
needed.  On output the diagonal and lower triangular part of the input
matrix <var>A</var> contain the matrix L. The upper triangular part
of the input matrix contains L^T, the diagonal terms being
identical for both L and L^T.  If the matrix is not
positive-definite then the decomposition will fail, returning the
error code <code>GSL_EDOM</code>.

        <p>When testing whether a matrix is positive-definite, disable the error
handler first to avoid triggering an error. 
</p></blockquote></div>

<div class="defun">
&mdash; Function: int <b>gsl_linalg_cholesky_solve</b> (<var>const gsl_matrix * cholesky, const gsl_vector * b, gsl_vector * x</var>)<var><a name="index-gsl_005flinalg_005fcholesky_005fsolve-1284"></a></var><br>
&mdash; Function: int <b>gsl_linalg_complex_cholesky_solve</b> (<var>const gsl_matrix_complex * cholesky, const gsl_vector_complex * b, gsl_vector_complex * x</var>)<var><a name="index-gsl_005flinalg_005fcomplex_005fcholesky_005fsolve-1285"></a></var><br>
<blockquote><p>These functions solve the system A x = b using the Cholesky
decomposition of A into the matrix <var>cholesky</var> given by
<code>gsl_linalg_cholesky_decomp</code> or
<code>gsl_linalg_complex_cholesky_decomp</code>. 
</p></blockquote></div>

<div class="defun">
&mdash; Function: int <b>gsl_linalg_cholesky_svx</b> (<var>const gsl_matrix * cholesky, gsl_vector * x</var>)<var><a name="index-gsl_005flinalg_005fcholesky_005fsvx-1286"></a></var><br>
&mdash; Function: int <b>gsl_linalg_complex_cholesky_svx</b> (<var>const gsl_matrix_complex * cholesky, gsl_vector_complex * x</var>)<var><a name="index-gsl_005flinalg_005fcomplex_005fcholesky_005fsvx-1287"></a></var><br>
<blockquote><p>These functions solve the system A x = b in-place using the
Cholesky decomposition of A into the matrix <var>cholesky</var> given
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 b, which is replaced by the solution on
output. 
</p></blockquote></div>

<hr>The GNU Scientific Library - a free numerical library licensed under the GNU GPL<br>Back to the <a href="/software/gsl/">GNU Scientific Library Homepage</a></body></html>