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
|
<!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 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 no
Invariant Sections and no cover texts. A copy of the license is
included in the section entitled "GNU Free Documentation License". -->
<!-- Created by GNU Texinfo 5.1, http://www.gnu.org/software/texinfo/ -->
<head>
<title>GNU Scientific Library – Reference Manual: LU Decomposition</title>
<meta name="description" content="GNU Scientific Library – Reference Manual: LU Decomposition">
<meta name="keywords" content="GNU Scientific Library – Reference Manual: LU 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="QR-Decomposition.html#QR-Decomposition" rel="next" title="QR Decomposition">
<link href="Linear-Algebra.html#Linear-Algebra" rel="previous" title="Linear Algebra">
<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="LU-Decomposition"></a>
<div class="header">
<p>
Next: <a href="QR-Decomposition.html#QR-Decomposition" accesskey="n" rel="next">QR 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="LU-Decomposition-1"></a>
<h3 class="section">14.1 LU Decomposition</h3>
<a name="index-LU-decomposition"></a>
<p>A general square matrix <em>A</em> has an <em>LU</em> decomposition into
upper and lower triangular matrices,
where <em>P</em> is a permutation matrix, <em>L</em> is unit lower
triangular matrix and <em>U</em> is upper triangular matrix. For square
matrices this decomposition can be used to convert the linear system
<em>A x = b</em> into a pair of triangular systems (<em>L y = P b</em>,
<em>U x = y</em>), which can be solved by forward and back-substitution.
Note that the <em>LU</em> decomposition is valid for singular matrices.
</p>
<dl>
<dt><a name="index-gsl_005flinalg_005fLU_005fdecomp"></a>Function: <em>int</em> <strong>gsl_linalg_LU_decomp</strong> <em>(gsl_matrix * <var>A</var>, gsl_permutation * <var>p</var>, int * <var>signum</var>)</em></dt>
<dt><a name="index-gsl_005flinalg_005fcomplex_005fLU_005fdecomp"></a>Function: <em>int</em> <strong>gsl_linalg_complex_LU_decomp</strong> <em>(gsl_matrix_complex * <var>A</var>, gsl_permutation * <var>p</var>, int * <var>signum</var>)</em></dt>
<dd><p>These functions factorize the square matrix <var>A</var> into the <em>LU</em>
decomposition <em>PA = LU</em>. On output the diagonal and upper
triangular part of the input matrix <var>A</var> contain the matrix
<em>U</em>. The lower triangular part of the input matrix (excluding the
diagonal) contains <em>L</em>. The diagonal elements of <em>L</em> are
unity, and are not stored.
</p>
<p>The permutation matrix <em>P</em> is encoded in the permutation
<var>p</var>. The <em>j</em>-th column of the matrix <em>P</em> is given by the
<em>k</em>-th column of the identity matrix, where <em>k = p_j</em> the
<em>j</em>-th element of the permutation vector. The sign of the
permutation is given by <var>signum</var>. It has the value <em>(-1)^n</em>,
where <em>n</em> is the number of interchanges in the permutation.
</p>
<p>The algorithm used in the decomposition is Gaussian Elimination with
partial pivoting (Golub & Van Loan, <cite>Matrix Computations</cite>,
Algorithm 3.4.1).
</p></dd></dl>
<a name="index-linear-systems_002c-solution-of"></a>
<dl>
<dt><a name="index-gsl_005flinalg_005fLU_005fsolve"></a>Function: <em>int</em> <strong>gsl_linalg_LU_solve</strong> <em>(const gsl_matrix * <var>LU</var>, const gsl_permutation * <var>p</var>, const gsl_vector * <var>b</var>, gsl_vector * <var>x</var>)</em></dt>
<dt><a name="index-gsl_005flinalg_005fcomplex_005fLU_005fsolve"></a>Function: <em>int</em> <strong>gsl_linalg_complex_LU_solve</strong> <em>(const gsl_matrix_complex * <var>LU</var>, const gsl_permutation * <var>p</var>, const gsl_vector_complex * <var>b</var>, gsl_vector_complex * <var>x</var>)</em></dt>
<dd><p>These functions solve the square system <em>A x = b</em> using the <em>LU</em>
decomposition of <em>A</em> into (<var>LU</var>, <var>p</var>) given by
<code>gsl_linalg_LU_decomp</code> or <code>gsl_linalg_complex_LU_decomp</code> as input.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005flinalg_005fLU_005fsvx"></a>Function: <em>int</em> <strong>gsl_linalg_LU_svx</strong> <em>(const gsl_matrix * <var>LU</var>, const gsl_permutation * <var>p</var>, gsl_vector * <var>x</var>)</em></dt>
<dt><a name="index-gsl_005flinalg_005fcomplex_005fLU_005fsvx"></a>Function: <em>int</em> <strong>gsl_linalg_complex_LU_svx</strong> <em>(const gsl_matrix_complex * <var>LU</var>, const gsl_permutation * <var>p</var>, gsl_vector_complex * <var>x</var>)</em></dt>
<dd><p>These functions solve the square system <em>A x = b</em> in-place using the
precomputed <em>LU</em> decomposition of <em>A</em> into (<var>LU</var>,<var>p</var>). 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>
<a name="index-refinement-of-solutions-in-linear-systems"></a>
<a name="index-iterative-refinement-of-solutions-in-linear-systems"></a>
<a name="index-linear-systems_002c-refinement-of-solutions"></a>
<dl>
<dt><a name="index-gsl_005flinalg_005fLU_005frefine"></a>Function: <em>int</em> <strong>gsl_linalg_LU_refine</strong> <em>(const gsl_matrix * <var>A</var>, const gsl_matrix * <var>LU</var>, const gsl_permutation * <var>p</var>, const gsl_vector * <var>b</var>, gsl_vector * <var>x</var>, gsl_vector * <var>residual</var>)</em></dt>
<dt><a name="index-gsl_005flinalg_005fcomplex_005fLU_005frefine"></a>Function: <em>int</em> <strong>gsl_linalg_complex_LU_refine</strong> <em>(const gsl_matrix_complex * <var>A</var>, const gsl_matrix_complex * <var>LU</var>, const gsl_permutation * <var>p</var>, const gsl_vector_complex * <var>b</var>, gsl_vector_complex * <var>x</var>, gsl_vector_complex * <var>residual</var>)</em></dt>
<dd><p>These functions apply an iterative improvement to <var>x</var>, the solution
of <em>A x = b</em>, from the precomputed <em>LU</em> decomposition of <em>A</em> into
(<var>LU</var>,<var>p</var>). The initial residual <em>r = A x - b</em> is also
computed and stored in <var>residual</var>.
</p></dd></dl>
<a name="index-inverse-of-a-matrix_002c-by-LU-decomposition"></a>
<a name="index-matrix-inverse"></a>
<dl>
<dt><a name="index-gsl_005flinalg_005fLU_005finvert"></a>Function: <em>int</em> <strong>gsl_linalg_LU_invert</strong> <em>(const gsl_matrix * <var>LU</var>, const gsl_permutation * <var>p</var>, gsl_matrix * <var>inverse</var>)</em></dt>
<dt><a name="index-gsl_005flinalg_005fcomplex_005fLU_005finvert"></a>Function: <em>int</em> <strong>gsl_linalg_complex_LU_invert</strong> <em>(const gsl_matrix_complex * <var>LU</var>, const gsl_permutation * <var>p</var>, gsl_matrix_complex * <var>inverse</var>)</em></dt>
<dd><p>These functions compute the inverse of a matrix <em>A</em> from its
<em>LU</em> decomposition (<var>LU</var>,<var>p</var>), storing the result in the
matrix <var>inverse</var>. The inverse is computed by solving the system
<em>A x = b</em> for each column of the identity matrix. It is preferable
to avoid direct use of the inverse whenever possible, as the linear
solver functions can obtain the same result more efficiently and
reliably (consult any introductory textbook on numerical linear algebra
for details).
</p></dd></dl>
<a name="index-determinant-of-a-matrix_002c-by-LU-decomposition"></a>
<a name="index-matrix-determinant"></a>
<dl>
<dt><a name="index-gsl_005flinalg_005fLU_005fdet"></a>Function: <em>double</em> <strong>gsl_linalg_LU_det</strong> <em>(gsl_matrix * <var>LU</var>, int <var>signum</var>)</em></dt>
<dt><a name="index-gsl_005flinalg_005fcomplex_005fLU_005fdet"></a>Function: <em>gsl_complex</em> <strong>gsl_linalg_complex_LU_det</strong> <em>(gsl_matrix_complex * <var>LU</var>, int <var>signum</var>)</em></dt>
<dd><p>These functions compute the determinant of a matrix <em>A</em> from its
<em>LU</em> decomposition, <var>LU</var>. The determinant is computed as the
product of the diagonal elements of <em>U</em> and the sign of the row
permutation <var>signum</var>.
</p></dd></dl>
<a name="index-logarithm-of-the-determinant-of-a-matrix"></a>
<dl>
<dt><a name="index-gsl_005flinalg_005fLU_005flndet"></a>Function: <em>double</em> <strong>gsl_linalg_LU_lndet</strong> <em>(gsl_matrix * <var>LU</var>)</em></dt>
<dt><a name="index-gsl_005flinalg_005fcomplex_005fLU_005flndet"></a>Function: <em>double</em> <strong>gsl_linalg_complex_LU_lndet</strong> <em>(gsl_matrix_complex * <var>LU</var>)</em></dt>
<dd><p>These functions compute the logarithm of the absolute value of the
determinant of a matrix <em>A</em>, <em>\ln|\det(A)|</em>, from its <em>LU</em>
decomposition, <var>LU</var>. This function may be useful if the direct
computation of the determinant would overflow or underflow.
</p></dd></dl>
<a name="index-sign-of-the-determinant-of-a-matrix"></a>
<dl>
<dt><a name="index-gsl_005flinalg_005fLU_005fsgndet"></a>Function: <em>int</em> <strong>gsl_linalg_LU_sgndet</strong> <em>(gsl_matrix * <var>LU</var>, int <var>signum</var>)</em></dt>
<dt><a name="index-gsl_005flinalg_005fcomplex_005fLU_005fsgndet"></a>Function: <em>gsl_complex</em> <strong>gsl_linalg_complex_LU_sgndet</strong> <em>(gsl_matrix_complex * <var>LU</var>, int <var>signum</var>)</em></dt>
<dd><p>These functions compute the sign or phase factor of the determinant of a
matrix <em>A</em>, <em>\det(A)/|\det(A)|</em>, from its <em>LU</em> decomposition,
<var>LU</var>.
</p></dd></dl>
<hr>
<div class="header">
<p>
Next: <a href="QR-Decomposition.html#QR-Decomposition" accesskey="n" rel="next">QR 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>
|