📄 linalg.texi
字号:
@cindex linear algebra@cindex solution of linear systems, Ax=b@cindex matrix factorization@cindex factorization of matricesThis chapter describes functions for solving linear systems. Thelibrary provides simple linear algebra operations which operate directlyon the @code{gsl_vector} and @code{gsl_matrix} objects. These areintended for use with ``small'' systems where simple algorithms areacceptable.@cindex LAPACK, recommended for linear algebraAnyone interested in large systems will want to use the sophisticatedroutines found in @sc{lapack}. The Fortran version of @sc{lapack} isrecommended as the standard package for large-scale linear algebra. Itsupports blocked algorithms, specialized data representations and otheroptimizations.The functions described in this chapter are declared in the header file@file{gsl_linalg.h}.@menu* LU Decomposition:: * QR Decomposition:: * QR Decomposition with Column Pivoting:: * Singular Value Decomposition:: * Cholesky Decomposition:: * Tridiagonal Decomposition of Real Symmetric Matrices:: * Tridiagonal Decomposition of Hermitian Matrices:: * Bidiagonalization:: * Householder Transformations:: * Householder solver for linear systems:: * Tridiagonal Systems:: * Linear Algebra Examples:: * Linear Algebra References and Further Reading:: @end menu@node LU Decomposition@section LU Decomposition@cindex LU decompositionA general square matrix @math{A} has an @math{LU} decomposition intoupper and lower triangular matrices,@tex\beforedisplay$$P A = L U$$\afterdisplay@end tex@ifinfo@exampleP A = L U@end example@end ifinfo@noindentwhere @math{P} is a permutation matrix, @math{L} is unit lowertriangular matrix and @math{U} is upper triangular matrix. For squarematrices this decomposition can be used to convert the linear system@math{A x = b} into a pair of triangular systems (@math{L y = P b},@math{U x = y}), which can be solved by forward and back-substitution.Note that the @math{LU} decomposition is valid for singular matrices.@deftypefun int gsl_linalg_LU_decomp (gsl_matrix * @var{A}, gsl_permutation * @var{p}, int * @var{signum})@deftypefunx int gsl_linalg_complex_LU_decomp (gsl_matrix_complex * @var{A}, gsl_permutation * @var{p}, int * @var{signum})These functions factorize the square matrix @var{A} into the @math{LU}decomposition @math{PA = LU}. On output the diagonal and uppertriangular part of the input matrix @var{A} contain the matrix@math{U}. The lower triangular part of the input matrix (excluding thediagonal) contains @math{L}. The diagonal elements of @math{L} areunity, and are not stored.The permutation matrix @math{P} is encoded in the permutation@var{p}. The @math{j}-th column of the matrix @math{P} is given by the@math{k}-th column of the identity matrix, where @math{k = p_j} the@math{j}-th element of the permutation vector. The sign of thepermutation is given by @var{signum}. It has the value @math{(-1)^n},where @math{n} is the number of interchanges in the permutation.The algorithm used in the decomposition is Gaussian Elimination withpartial pivoting (Golub & Van Loan, @cite{Matrix Computations},Algorithm 3.4.1).@end deftypefun@cindex linear systems, solution of@deftypefun int gsl_linalg_LU_solve (const gsl_matrix * @var{LU}, const gsl_permutation * @var{p}, const gsl_vector * @var{b}, gsl_vector * @var{x})@deftypefunx int gsl_linalg_complex_LU_solve (const gsl_matrix_complex * @var{LU}, const gsl_permutation * @var{p}, const gsl_vector_complex * @var{b}, gsl_vector_complex * @var{x})These functions solve the square system @math{A x = b} using the @math{LU}decomposition of @math{A} into (@var{LU}, @var{p}) given by@code{gsl_linalg_LU_decomp} or @code{gsl_linalg_complex_LU_decomp}.@end deftypefun@deftypefun int gsl_linalg_LU_svx (const gsl_matrix * @var{LU}, const gsl_permutation * @var{p}, gsl_vector * @var{x})@deftypefunx int gsl_linalg_complex_LU_svx (const gsl_matrix_complex * @var{LU}, const gsl_permutation * @var{p}, gsl_vector_complex * @var{x})These functions solve the square system @math{A x = b} in-place using the@math{LU} decomposition of @math{A} into (@var{LU},@var{p}). On input@var{x} should contain the right-hand side @math{b}, which is replacedby the solution on output.@end deftypefun@cindex refinement of solutions in linear systems@cindex iterative refinement of solutions in linear systems@cindex linear systems, refinement of solutions@deftypefun int gsl_linalg_LU_refine (const gsl_matrix * @var{A}, const gsl_matrix * @var{LU}, const gsl_permutation * @var{p}, const gsl_vector * @var{b}, gsl_vector * @var{x}, gsl_vector * @var{residual})@deftypefunx int gsl_linalg_complex_LU_refine (const gsl_matrix_complex * @var{A}, const gsl_matrix_complex * @var{LU}, const gsl_permutation * @var{p}, const gsl_vector_complex * @var{b}, gsl_vector_complex * @var{x}, gsl_vector_complex * @var{residual})These functions apply an iterative improvement to @var{x}, the solutionof @math{A x = b}, using the @math{LU} decomposition of @math{A} into(@var{LU},@var{p}). The initial residual @math{r = A x - b} is alsocomputed and stored in @var{residual}.@end deftypefun@cindex inverse of a matrix, by LU decomposition@cindex matrix inverse@deftypefun int gsl_linalg_LU_invert (const gsl_matrix * @var{LU}, const gsl_permutation * @var{p}, gsl_matrix * @var{inverse})@deftypefunx int gsl_linalg_complex_LU_invert (const gsl_matrix_complex * @var{LU}, const gsl_permutation * @var{p}, gsl_matrix_complex * @var{inverse})These functions compute the inverse of a matrix @math{A} from its@math{LU} decomposition (@var{LU},@var{p}), storing the result in thematrix @var{inverse}. The inverse is computed by solving the system@math{A x = b} for each column of the identity matrix. It is preferableto avoid direct use of the inverse whenever possible, as the linearsolver functions can obtain the same result more efficiently andreliably (consult any introductory textbook on numerical linear algebrafor details).@end deftypefun@cindex determinant of a matrix, by LU decomposition@cindex matrix determinant@deftypefun double gsl_linalg_LU_det (gsl_matrix * @var{LU}, int @var{signum})@deftypefunx gsl_complex gsl_linalg_complex_LU_det (gsl_matrix_complex * @var{LU}, int @var{signum})These functions compute the determinant of a matrix @math{A} from its@math{LU} decomposition, @var{LU}. The determinant is computed as theproduct of the diagonal elements of @math{U} and the sign of the rowpermutation @var{signum}.@end deftypefun@cindex logarithm of the determinant of a matrix@deftypefun double gsl_linalg_LU_lndet (gsl_matrix * @var{LU})@deftypefunx double gsl_linalg_complex_LU_lndet (gsl_matrix_complex * @var{LU})These functions compute the logarithm of the absolute value of thedeterminant of a matrix @math{A}, @math{\ln|\det(A)|}, from its @math{LU}decomposition, @var{LU}. This function may be useful if the directcomputation of the determinant would overflow or underflow.@end deftypefun@cindex sign of the determinant of a matrix@deftypefun int gsl_linalg_LU_sgndet (gsl_matrix * @var{LU}, int @var{signum})@deftypefunx gsl_complex gsl_linalg_complex_LU_sgndet (gsl_matrix_complex * @var{LU}, int @var{signum})These functions compute the sign or phase factor of the determinant of amatrix @math{A}, @math{\det(A)/|\det(A)|}, from its @math{LU} decomposition,@var{LU}.@end deftypefun@node QR Decomposition@section QR Decomposition@cindex QR decompositionA general rectangular @math{M}-by-@math{N} matrix @math{A} has a@math{QR} decomposition into the product of an orthogonal@math{M}-by-@math{M} square matrix @math{Q} (where @math{Q^T Q = I}) andan @math{M}-by-@math{N} right-triangular matrix @math{R},@tex\beforedisplay$$A = Q R$$\afterdisplay@end tex@ifinfo@exampleA = Q R@end example@end ifinfo@noindentThis decomposition can be used to convert the linear system @math{A x =b} into the triangular system @math{R x = Q^T b}, which can be solved byback-substitution. Another use of the @math{QR} decomposition is tocompute an orthonormal basis for a set of vectors. The first @math{N}columns of @math{Q} form an orthonormal basis for the range of @math{A},@math{ran(A)}, when @math{A} has full column rank.@deftypefun int gsl_linalg_QR_decomp (gsl_matrix * @var{A}, gsl_vector * @var{tau})This function factorizes the @math{M}-by-@math{N} matrix @var{A} intothe @math{QR} decomposition @math{A = Q R}. On output the diagonal andupper triangular part of the input matrix contain the matrix@math{R}. The vector @var{tau} and the columns of the lower triangularpart of the matrix @var{A} contain the Householder coefficients andHouseholder vectors which encode the orthogonal matrix @var{Q}. Thevector @var{tau} must be of length @math{k=\min(M,N)}. The matrix@math{Q} is related to these components by, @math{Q = Q_k ... Q_2 Q_1}where @math{Q_i = I - \tau_i v_i v_i^T} and @math{v_i} is theHouseholder vector @math{v_i =(0,...,1,A(i+1,i),A(i+2,i),...,A(m,i))}. This is the same storage schemeas used by @sc{lapack}.The algorithm used to perform the decomposition is Householder QR (Golub& Van Loan, @cite{Matrix Computations}, Algorithm 5.2.1).@end deftypefun@deftypefun int gsl_linalg_QR_solve (const gsl_matrix * @var{QR}, const gsl_vector * @var{tau}, const gsl_vector * @var{b}, gsl_vector * @var{x})This function solves the square system @math{A x = b} using the @math{QR}decomposition of @math{A} into (@var{QR}, @var{tau}) given by@code{gsl_linalg_QR_decomp}. The least-squares solution for rectangular systems canbe found using @code{gsl_linalg_QR_lssolve}.@end deftypefun@deftypefun int gsl_linalg_QR_svx (const gsl_matrix * @var{QR}, const gsl_vector * @var{tau}, gsl_vector * @var{x})This function solves the square system @math{A x = b} in-place using the@math{QR} decomposition of @math{A} into (@var{QR},@var{tau}) given by@code{gsl_linalg_QR_decomp}. On input @var{x} should contain theright-hand side @math{b}, which is replaced by the solution on output.@end deftypefun@deftypefun int gsl_linalg_QR_lssolve (const gsl_matrix * @var{QR}, const gsl_vector * @var{tau}, const gsl_vector * @var{b}, gsl_vector * @var{x}, gsl_vector * @var{residual})This function finds the least squares solution to the overdeterminedsystem @math{A x = b} where the matrix @var{A} has more rows thancolumns. The least squares solution minimizes the Euclidean norm of theresidual, @math{||Ax - b||}.The routine uses the @math{QR} decompositionof @math{A} into (@var{QR}, @var{tau}) given by@code{gsl_linalg_QR_decomp}. The solution is returned in @var{x}. Theresidual is computed as a by-product and stored in @var{residual}.@end deftypefun@deftypefun int gsl_linalg_QR_QTvec (const gsl_matrix * @var{QR}, const gsl_vector * @var{tau}, gsl_vector * @var{v})This function applies the matrix @math{Q^T} encoded in the decomposition(@var{QR},@var{tau}) to the vector @var{v}, storing the result @math{Q^Tv} in @var{v}. The matrix multiplication is carried out directly usingthe encoding of the Householder vectors without needing to form the fullmatrix @math{Q^T}.@end deftypefun@deftypefun int gsl_linalg_QR_Qvec (const gsl_matrix * @var{QR}, const gsl_vector * @var{tau}, gsl_vector * @var{v})This function applies the matrix @math{Q} encoded in the decomposition(@var{QR},@var{tau}) to the vector @var{v}, storing the result @math{Qv} in @var{v}. The matrix multiplication is carried out directly usingthe encoding of the Householder vectors without needing to form the fullmatrix @math{Q}.@end deftypefun@deftypefun int gsl_linalg_QR_Rsolve (const gsl_matrix * @var{QR}, const gsl_vector * @var{b}, gsl_vector * @var{x})This function solves the triangular system @math{R x = b} for@var{x}. It may be useful if the product @math{b' = Q^T b} has alreadybeen computed using @code{gsl_linalg_QR_QTvec}.@end deftypefun@deftypefun int gsl_linalg_QR_Rsvx (const gsl_matrix * @var{QR}, gsl_vector * @var{x})This function solves the triangular system @math{R x = b} for @var{x}in-place. On input @var{x} should contain the right-hand side @math{b}and is replaced by the solution on output. This function may be useful ifthe product @math{b' = Q^T b} has already been computed using@code{gsl_linalg_QR_QTvec}.@end deftypefun@deftypefun int gsl_linalg_QR_unpack (const gsl_matrix * @var{QR}, const gsl_vector * @var{tau}, gsl_matrix * @var{Q}, gsl_matrix * @var{R})This function unpacks the encoded @math{QR} decomposition(@var{QR},@var{tau}) into the matrices @var{Q} and @var{R}, where@var{Q} is @math{M}-by-@math{M} and @var{R} is @math{M}-by-@math{N}. @end deftypefun@deftypefun int gsl_linalg_QR_QRsolve (gsl_matrix * @var{Q}, gsl_matrix * @var{R}, const gsl_vector * @var{b}, gsl_vector * @var{x})This function solves the system @math{R x = Q^T b} for @var{x}. It canbe used when the @math{QR} decomposition of a matrix is available inunpacked form as (@var{Q}, @var{R}).@end deftypefun@deftypefun int gsl_linalg_QR_update (gsl_matrix * @var{Q}, gsl_matrix * @var{R}, gsl_vector * @var{w}, const gsl_vector * @var{v})This function performs a rank-1 update @math{w v^T} of the @math{QR}decomposition (@var{Q}, @var{R}). The update is given by @math{Q'R' = QR + w v^T} where the output matrices @math{Q'} and @math{R'} are alsoorthogonal and right triangular. Note that @var{w} is destroyed by theupdate.@end deftypefun@deftypefun int gsl_linalg_R_solve (const gsl_matrix * @var{R}, const gsl_vector * @var{b}, gsl_vector * @var{x})This function solves the triangular system @math{R x = b} for the@math{N}-by-@math{N} matrix @var{R}.@end deftypefun@deftypefun int gsl_linalg_R_svx (const gsl_matrix * @var{R}, gsl_vector * @var{x})This function solves the triangular system @math{R x = b} in-place. Oninput @var{x} should contain the right-hand side @math{b}, which isreplaced by the solution on output.@end deftypefun@node QR Decomposition with Column Pivoting@section QR Decomposition with Column Pivoting@cindex QR decomposition with column pivotingThe @math{QR} decomposition can be extended to the rank deficient caseby introducing a column permutation @math{P},@tex\beforedisplay$$A P = Q R$$\afterdisplay@end tex@ifinfo@exampleA P = Q R@end example@end ifinfo@noindentThe first @math{r} columns of @math{Q} form an orthonormal basisfor the range of @math{A} for a matrix with column rank @math{r}. Thisdecomposition can also be used to convert the linear system @math{A x =b} into the triangular system @math{R y = Q^T b, x = P y}, which can besolved by back-substitution and permutation. We denote the @math{QR}decomposition with column pivoting by @math{QRP^T} since @math{A = Q RP^T}.@deftypefun int gsl_linalg_QRPT_decomp (gsl_matrix * @var{A}, gsl_vector * @var{tau}, gsl_permutation * @var{p}, int * @var{signum}, gsl_vector * @var{norm})This function factorizes the @math{M}-by-@math{N} matrix @var{A} intothe @math{QRP^T} decomposition @math{A = Q R P^T}. On output thediagonal and upper triangular part of the input matrix contain thematrix @math{R}. The permutation matrix @math{P} is stored in thepermutation @var{p}. The sign of the permutation is given by@var{signum}. It has the value @math{(-1)^n}, where @math{n} is thenumber of interchanges in the permutation. The vector @var{tau} and thecolumns of the lower triangular part of the matrix @var{A} contain theHouseholder coefficients and vectors which encode the orthogonal matrix@var{Q}. The vector @var{tau} must be of length @math{k=\min(M,N)}. Thematrix @math{Q} is related to these components by, @math{Q = Q_k ... Q_2Q_1} where @math{Q_i = I - \tau_i v_i v_i^T} and @math{v_i} is theHouseholder vector @math{v_i =(0,...,1,A(i+1,i),A(i+2,i),...,A(m,i))}. This is the same storage schemeas used by @sc{lapack}. The vector @var{norm} is a workspace of length@var{N} used for column pivoting.The algorithm used to perform the decomposition is Householder QR withcolumn pivoting (Golub & Van Loan, @cite{Matrix Computations}, Algorithm5.4.1).@end deftypefun
⌨️ 快捷键说明
复制代码
Ctrl + C
搜索代码
Ctrl + F
全屏模式
F11
切换主题
Ctrl + Shift + D
显示快捷键
?
增大字号
Ctrl + =
减小字号
Ctrl + -