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
|
<html lang="en">
<head>
<title>Tridiagonal Decomposition of Real Symmetric Matrices - 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="Cholesky-Decomposition.html" title="Cholesky Decomposition">
<link rel="next" href="Tridiagonal-Decomposition-of-Hermitian-Matrices.html" title="Tridiagonal Decomposition of Hermitian 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="Tridiagonal-Decomposition-of-Real-Symmetric-Matrices"></a>
Next: <a rel="next" accesskey="n" href="Tridiagonal-Decomposition-of-Hermitian-Matrices.html">Tridiagonal Decomposition of Hermitian Matrices</a>,
Previous: <a rel="previous" accesskey="p" href="Cholesky-Decomposition.html">Cholesky Decomposition</a>,
Up: <a rel="up" accesskey="u" href="Linear-Algebra.html">Linear Algebra</a>
<hr>
</div>
<h3 class="section">13.6 Tridiagonal Decomposition of Real Symmetric Matrices</h3>
<p><a name="index-tridiagonal-decomposition-1288"></a>
A symmetric matrix A can be factorized by similarity
transformations into the form,
<pre class="example"> A = Q T Q^T
</pre>
<p class="noindent">where Q is an orthogonal matrix and T is a symmetric
tridiagonal matrix.
<div class="defun">
— Function: int <b>gsl_linalg_symmtd_decomp</b> (<var>gsl_matrix * A, gsl_vector * tau</var>)<var><a name="index-gsl_005flinalg_005fsymmtd_005fdecomp-1289"></a></var><br>
<blockquote><p>This function factorizes the symmetric square matrix <var>A</var> into the
symmetric tridiagonal decomposition Q T Q^T. On output the
diagonal and subdiagonal part of the input matrix <var>A</var> contain the
tridiagonal matrix T. The remaining lower triangular part of the
input matrix contains the Householder vectors which, together with the
Householder coefficients <var>tau</var>, encode the orthogonal matrix
Q. This storage scheme is the same as used by <span class="sc">lapack</span>. The
upper triangular part of <var>A</var> is not referenced.
</p></blockquote></div>
<div class="defun">
— Function: int <b>gsl_linalg_symmtd_unpack</b> (<var>const gsl_matrix * A, const gsl_vector * tau, gsl_matrix * Q, gsl_vector * diag, gsl_vector * subdiag</var>)<var><a name="index-gsl_005flinalg_005fsymmtd_005funpack-1290"></a></var><br>
<blockquote><p>This function unpacks the encoded symmetric tridiagonal decomposition
(<var>A</var>, <var>tau</var>) obtained from <code>gsl_linalg_symmtd_decomp</code> into
the orthogonal matrix <var>Q</var>, the vector of diagonal elements <var>diag</var>
and the vector of subdiagonal elements <var>subdiag</var>.
</p></blockquote></div>
<div class="defun">
— Function: int <b>gsl_linalg_symmtd_unpack_T</b> (<var>const gsl_matrix * A, gsl_vector * diag, gsl_vector * subdiag</var>)<var><a name="index-gsl_005flinalg_005fsymmtd_005funpack_005fT-1291"></a></var><br>
<blockquote><p>This function unpacks the diagonal and subdiagonal of the encoded
symmetric tridiagonal decomposition (<var>A</var>, <var>tau</var>) obtained from
<code>gsl_linalg_symmtd_decomp</code> into the vectors <var>diag</var> and <var>subdiag</var>.
</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>
|